source: trunk/LMDZ.MARS/libf/aeronomars/moldiffcoeff.F @ 1242

Last change on this file since 1242 was 1047, checked in by emillour, 11 years ago

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File size: 6.5 KB
Line 
1      subroutine moldiffcoeff(dij)
2
3      use tracer_mod, only: igcm_co2, igcm_co, igcm_o, igcm_o1d,
4     &                      igcm_o2, igcm_o3, igcm_h, igcm_h2, igcm_oh,
5     &                      igcm_ho2, igcm_h2o2, igcm_n2, igcm_ar,
6     &                      igcm_h2o_vap, mmol
7
8       IMPLICIT NONE
9c=======================================================================
10c   subject:
11c   --------
12c   Computing molecular diffusion coefficients
13c   following Nair 94 (pg 131)
14c   author:  MAC 2002
15c   ------
16c
17c=======================================================================
18#include "dimensions.h"
19#include "dimphys.h"
20#include "callkeys.h"
21!#include "comdiurn.h"
22#include "chimiedata.h"
23!#include "tracer.h"
24!#include "conc.h"
25
26c-----------------------------------------------------------------------
27c    Input/Output
28c    ------------
29       integer,parameter :: ncompmoldiff = 14
30      real dij(ncompmoldiff,ncompmoldiff)
31
32c    Local variables:
33c    ---------------
34      INTEGER nq, n, nn, i,iq
35cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
36c     tracer numbering in the molecular diffusion
37cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
38c  Atomic oxygen must always be the LAST species of the list as
39c it is the dominant species at high altitudes. 
40      integer,parameter :: i_co   = 1
41      integer,parameter :: i_n2   = 2
42      integer,parameter :: i_o2   = 3
43      integer,parameter :: i_co2  = 4
44      integer,parameter :: i_h2   = 5
45      integer,parameter :: i_h    = 6
46      integer,parameter :: i_oh   = 7
47      integer,parameter :: i_ho2  = 8
48      integer,parameter :: i_h2o  = 9
49      integer,parameter :: i_h2o2 = 10
50      integer,parameter :: i_o1d  = 11
51      integer,parameter :: i_o3   = 12
52      integer,parameter :: i_ar   = 13
53      integer,parameter :: i_o    = 14
54
55! Tracer indexes in the GCM:
56      integer,save :: g_co2=0
57      integer,save :: g_co=0
58      integer,save :: g_o=0
59      integer,save :: g_o1d=0
60      integer,save :: g_o2=0
61      integer,save :: g_o3=0
62      integer,save :: g_h=0
63      integer,save :: g_h2=0
64      integer,save :: g_oh=0
65      integer,save :: g_ho2=0
66      integer,save :: g_h2o2=0
67      integer,save :: g_n2=0
68      integer,save :: g_ar=0
69      integer,save :: g_h2o=0
70
71      integer,save :: gcmind(ncompmoldiff)
72
73      real dnh
74      logical,save :: firstcall=.true.
75
76! Initializations at first call (and some sanity checks)
77      if (firstcall) then
78        ! identify the indexes of the tracers we'll need
79        g_co2=igcm_co2
80        if (g_co2.eq.0) then
81          write(*,*) "moldiffcoeff: Error; no CO2 tracer !!!"
82          stop
83        endif
84        g_co=igcm_co
85        if (g_co.eq.0) then
86          write(*,*) "moldiffcoeff: Error; no CO tracer !!!"
87          stop
88        endif
89        g_o=igcm_o
90        if (g_o.eq.0) then
91          write(*,*) "moldiffcoeff: Error; no O tracer !!!"
92          stop
93        endif
94        g_o1d=igcm_o1d
95        if (g_o1d.eq.0) then
96          write(*,*) "moldiffcoeff: Error; no O1D tracer !!!"
97          stop
98        endif
99        g_o2=igcm_o2
100        if (g_o2.eq.0) then
101          write(*,*) "moldiffcoeff: Error; no O2 tracer !!!"
102          stop
103        endif
104        g_o3=igcm_o3
105        if (g_o3.eq.0) then
106          write(*,*) "moldiffcoeff: Error; no O3 tracer !!!"
107          stop
108        endif
109        g_h=igcm_h
110        if (g_h.eq.0) then
111          write(*,*) "moldiffcoeff: Error; no H tracer !!!"
112          stop
113        endif
114        g_h2=igcm_h2
115        if (g_h2.eq.0) then
116          write(*,*) "moldiffcoeff: Error; no H2 tracer !!!"
117          stop
118        endif
119        g_oh=igcm_oh
120        if (g_oh.eq.0) then
121          write(*,*) "moldiffcoeff: Error; no OH tracer !!!"
122          stop
123        endif
124        g_ho2=igcm_ho2
125        if (g_ho2.eq.0) then
126          write(*,*) "moldiffcoeff: Error; no HO2 tracer !!!"
127          stop
128        endif
129        g_h2o2=igcm_h2o2
130        if (g_h2o2.eq.0) then
131          write(*,*) "moldiffcoeff: Error; no H2O2 tracer !!!"
132          stop
133        endif
134        g_n2=igcm_n2
135        if (g_n2.eq.0) then
136          write(*,*) "moldiffcoeff: Error; no N2 tracer !!!"
137          stop
138        endif
139        g_ar=igcm_ar
140        if (g_ar.eq.0) then
141          write(*,*) "moldiffcoeff: Error; no AR tracer !!!"
142          stop
143        endif
144        g_h2o=igcm_h2o_vap
145        if (g_h2o.eq.0) then
146          write(*,*) "moldiffcoeff: Error; no water vapor tracer !!!"
147          stop
148        endif
149
150c
151cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
152c    fill array to relate local indexes to gcm indexes
153cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
154
155        gcmind(i_co)  =   g_co
156        gcmind(i_n2)  =   g_n2
157        gcmind(i_o2)  =   g_o2
158        gcmind(i_co2) =   g_co2
159        gcmind(i_h2)  =   g_h2
160        gcmind(i_h)   =   g_h
161        gcmind(i_oh)  =   g_oh
162        gcmind(i_ho2) =   g_ho2
163        gcmind(i_h2o) =   g_h2o
164        gcmind(i_h2o2)=   g_h2o2
165        gcmind(i_o1d) =   g_o1d
166        gcmind(i_o3)  =   g_o3
167        gcmind(i_o)   =   g_o
168        gcmind(i_ar)   =  g_ar
169c
170cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
171        firstcall= .false.
172      endif ! of if (firstcall)
173
174
175      dij(i_h2,i_co)   = 0.0000651
176      dij(i_h2,i_n2)   = 0.0000674
177      dij(i_h2,i_o2)   = 0.0000697
178      dij(i_h2,i_co2)  = 0.0000550
179      dij(i_h2,i_h2)   = 0.0
180      dij(i_h2,i_h)    = 0.0
181      dij(i_h2,i_oh)   = 0.0    !0003
182      dij(i_h2,i_ho2)  = 0.0    !0003
183      dij(i_h2,i_h2o)  = 0.0    !0003
184      dij(i_h2,i_h2o2) = 0.0    !0003
185      dij(i_h2,i_o1d)  = 0.0
186      dij(i_h2,i_o3)   = 0.0    !0003
187      dij(i_h2,i_o)    = 0.0
188      dij(i_h2,i_ar)   = 0.0
189
190c      dij(i_h,i_o)     = 0.0000144
191      dij(i_h,i_o)     = 0.000114
192
193       print*,'moldiffcoef: COEFF CALC'
194       open(56,file='coeffs.dat',status='unknown')
195      do n=1,ncompmoldiff
196        if (dij(i_h2,n).gt.0.0) then
197          do nn=n,ncompmoldiff
198            dij(nn,n)=dij(i_h2,n)
199     &                  *sqrt(mmol(g_h2)/mmol(gcmind(nn)))
200            if(n.eq.nn) dij(nn,n)=1.0
201            dij(n,nn)=dij(nn,n)
202          enddo
203        endif
204        if (dij(i_h2,n).eq.0.0) then
205          dnh=dij(i_h,i_o)*sqrt(mmol(g_o)/mmol(gcmind(n)))
206          do nn=n,ncompmoldiff
207            dij(nn,n)=dnh*sqrt(mmol(g_h)/mmol(gcmind(nn)))
208            if(n.eq.nn) dij(nn,n)=1.0
209            dij(n,nn)=dij(nn,n)
210          enddo
211        endif
212      enddo
213
214      do n=1,ncompmoldiff
215        do nn=n,ncompmoldiff
216          write(56,*) n,nn,dij(n,nn)    !*1.e5/1.381e-23/(273**1.75)
217        enddo
218      enddo
219      close(56)
220
221
222      return   
223      end
Note: See TracBrowser for help on using the repository browser.