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

Last change on this file since 1621 was 1266, checked in by aslmd, 11 years ago

LMDZ.MARS
IMPORTANT CHANGE

  • Remove all reference/use of nlayermx and dimphys.h
  • Made use of automatic arrays whenever arrays are needed with dimension nlayer
  • Remove lots of obsolete reference to dimensions.h
  • Converted iono.h and param_v4.h into corresponding modules

(with embedded subroutine to allocate arrays)
(no arrays allocated if thermosphere not used)

  • Deleted param.h and put contents into module param_v4_h
  • Adapted testphys1d, newstart, etc...
  • Made DATA arrays in param_read to be initialized by subroutine

fill_data_thermos in module param_v4_h

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