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

Last change on this file since 2325 was 2030, checked in by flefevre, 6 years ago

photolyse on-line:

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