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

Last change on this file since 3807 was 3726, checked in by emillour, 2 months ago

Mars PCM:
Turn "callkeys.h" into module "callkeys_mod.F90"
EM

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