source: trunk/LMDZ.MARS/libf/aeronomars/concentrations.F @ 1198

Last change on this file since 1198 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: 8.7 KB
Line 
1      SUBROUTINE concentrations(ngrid,nlayer,nq,
2     &                          pplay,pt,pdt,pq,pdq,ptimestep)
3                                             
4      use tracer_mod, only: igcm_co2, igcm_co, igcm_o, igcm_o1d,
5     &                      igcm_o2, igcm_o3, igcm_h, igcm_h2,
6     &                      igcm_oh, igcm_ho2, igcm_n2, igcm_ar,
7     &                      igcm_h2o_vap, igcm_n, igcm_no, igcm_no2,
8     &                      igcm_n2d, igcm_co2plus, igcm_oplus,
9     &                      igcm_o2plus, igcm_coplus, igcm_cplus,
10     &                      igcm_nplus, igcm_noplus, igcm_n2plus,
11     &                      igcm_hplus, igcm_hco2plus, mmol
12      use conc_mod, only: mmean, Akknew, rnew, cpnew
13                         
14      implicit none
15
16!=======================================================================
17! CALCULATION OF MEAN MOLECULAR MASS, Cp, Akk and R
18!
19! mmean(ngrid,nlayer)   amu
20! cpnew(ngrid,nlayer)   J/kg/K
21! rnew(ngrid,nlayer)    J/kg/K
22! akknew(ngrid,nlayer)  coefficient of thermal concduction
23!
24! version: April 2012 - Franck Lefevre
25!=======================================================================
26
27!     declarations
28 
29!#include "dimensions.h"
30!#include "dimphys.h"
31#include "comcstfi.h"
32#include "callkeys.h"
33!#include "comdiurn.h"
34#include "chimiedata.h"
35!#include "tracer.h"
36!#include "conc.h"
37
38!     input/output
39
40      integer,intent(in) :: ngrid ! number of atmospheric columns
41      integer,intent(in) :: nlayer ! number of atmospheric layers
42      integer,intent(in) :: nq ! number of tracers
43      real,intent(in) :: pplay(ngrid,nlayer)
44      real,intent(in) :: pt(ngrid,nlayer)
45      real,intent(in) :: pdt(ngrid,nlayer)
46      real,intent(in) :: pq(ngrid,nlayer,nq)
47      real,intent(in) :: pdq(ngrid,nlayer,nq)
48      real,intent(in) :: ptimestep
49
50!     local variables
51
52      integer       :: i, l, ig, iq
53      integer, save :: nbq
54      integer,allocatable,save :: niq(:)
55      real          :: ni(nq), ntot
56      real          :: zq(ngrid, nlayer, nq)
57      real          :: zt(ngrid, nlayer)
58      real,allocatable,save    :: aki(:)
59      real,allocatable,save    :: cpi(:)
60
61      logical, save :: firstcall = .true.
62
63      if (firstcall) then
64
65         ! allocate local saved arrays:
66         allocate(aki(nq))
67         allocate(cpi(nq))
68         allocate(niq(nq))
69!        find index of chemical tracers to use
70!        initialize thermal conductivity and specific heat coefficients
71!        !? values are estimated
72
73         nbq = 0 ! to count number of tracers used in this subroutine
74
75         if (igcm_co2 /= 0) then
76            nbq = nbq + 1
77            niq(nbq) = igcm_co2
78            aki(nbq) = 3.072e-4
79            cpi(nbq) = 0.834e3
80         end if
81         if (igcm_co /= 0) then
82            nbq = nbq + 1
83            niq(nbq) = igcm_co
84            aki(nbq) = 4.87e-4
85            cpi(nbq) = 1.034e3
86         end if
87         if (igcm_o /= 0) then
88            nbq = nbq + 1
89            niq(nbq) = igcm_o
90            aki(nbq) = 7.59e-4
91            cpi(nbq) = 1.3e3
92         end if
93         if (igcm_o1d /= 0) then
94            nbq = nbq + 1
95            niq(nbq) = igcm_o1d
96            aki(nbq) = 7.59e-4  !?
97            cpi(nbq) = 1.3e3    !?
98         end if
99         if (igcm_o2 /= 0) then
100            nbq = nbq + 1
101            niq(nbq) = igcm_o2
102            aki(nbq) = 5.68e-4
103            cpi(nbq) = 0.9194e3
104         end if
105         if (igcm_o3 /= 0) then
106            nbq = nbq + 1
107            niq(nbq) = igcm_o3
108            aki(nbq) = 3.00e-4  !?
109            cpi(nbq) = 0.800e3  !?
110         end if
111         if (igcm_h /= 0) then
112            nbq = nbq + 1
113            niq(nbq) = igcm_h
114            aki(nbq) = 0.0
115            cpi(nbq) = 20.780e3
116         end if
117         if (igcm_h2 /= 0) then
118            nbq = nbq + 1
119            niq(nbq) = igcm_h2
120            aki(nbq) = 36.314e-4
121            cpi(nbq) = 14.266e3
122         end if
123         if (igcm_oh /= 0) then
124            nbq = nbq + 1
125            niq(nbq) = igcm_oh
126            aki(nbq)  = 7.00e-4 !?
127            cpi(nbq)  = 1.045e3
128         end if
129         if (igcm_ho2 /= 0) then
130            nbq = nbq + 1
131            niq(nbq) = igcm_ho2
132            aki(nbq) = 0.0
133            cpi(nbq) = 1.065e3  !?
134         end if
135         if (igcm_n2 /= 0) then
136            nbq = nbq + 1
137            niq(nbq) = igcm_n2
138            aki(nbq) = 5.6e-4
139            cpi(nbq) = 1.034e3
140         end if
141         if (igcm_ar /= 0) then
142            nbq = nbq + 1
143            niq(nbq) = igcm_ar
144            aki(nbq) = 0.0      !?
145            cpi(nbq) = 1.000e3  !?
146         end if
147         if (igcm_h2o_vap /= 0) then
148            nbq = nbq + 1
149            niq(nbq) = igcm_h2o_vap
150            aki(nbq) = 0.0
151            cpi(nbq) = 1.870e3
152         end if
153         if (igcm_n /= 0) then
154            nbq = nbq + 1
155            niq(nbq) = igcm_n
156            aki(nbq) = 0.0
157            cpi(nbq) = 0.0
158         endif
159         if(igcm_no /= 0) then
160            nbq = nbq + 1
161            niq(nbq) = igcm_no
162            aki(nbq) = 0.0
163            cpi(nbq) = 0.0
164         endif
165         if(igcm_no2 /= 0) then
166            nbq = nbq + 1
167            niq(nbq) = igcm_no2
168            aki(nbq) = 0.0
169            cpi(nbq) = 0.0
170         endif
171         if(igcm_n2d /= 0) then
172            nbq = nbq + 1
173            niq(nbq) = igcm_n2d
174            aki(nbq) = 0.0
175            cpi(nbq) = 0.0
176         endif
177         if(igcm_co2plus /= 0) then
178            nbq = nbq + 1
179            niq(nbq) = igcm_co2plus
180            aki(nbq) = 0.0
181            cpi(nbq) = 0.0
182         endif
183         if(igcm_oplus /= 0) then
184            nbq = nbq + 1
185            niq(nbq) = igcm_oplus
186            aki(nbq) = 0.0
187            cpi(nbq) = 0.0
188         endif
189         if(igcm_o2plus /= 0) then
190            nbq = nbq + 1
191            niq(nbq) = igcm_o2plus
192            aki(nbq) = 0.0
193            cpi(nbq) = 0.0
194         endif
195         if(igcm_coplus /= 0) then
196            nbq = nbq + 1
197            niq(nbq) = igcm_coplus
198            aki(nbq) = 0.0
199            cpi(nbq) = 0.0
200         endif
201         if(igcm_cplus /= 0) then
202            nbq = nbq + 1
203            niq(nbq) = igcm_cplus
204            aki(nbq) = 0.0
205            cpi(nbq) = 0.0
206         endif
207         if(igcm_nplus /= 0) then
208            nbq = nbq + 1
209            niq(nbq) = igcm_nplus
210            aki(nbq) = 0.0
211            cpi(nbq) = 0.0
212         endif
213         if(igcm_noplus /= 0) then
214            nbq = nbq + 1
215            niq(nbq) = igcm_noplus
216            aki(nbq) = 0.0
217            cpi(nbq) = 0.0
218         endif
219         if(igcm_n2plus /= 0) then
220            nbq = nbq + 1
221            niq(nbq) = igcm_n2plus
222            aki(nbq) = 0.0
223            cpi(nbq) = 0.0
224         endif
225         if(igcm_hplus /= 0) then
226            nbq = nbq + 1
227            niq(nbq) = igcm_hplus
228            aki(nbq) = 0.0
229            cpi(nbq) = 0.0
230         endif
231         if(igcm_hco2plus /= 0) then
232            nbq = nbq + 1
233            niq(nbq) = igcm_hco2plus
234            aki(nbq) = 0.0
235            cpi(nbq) = 0.0
236         endif
237         
238         ! tell the world about it:
239         write(*,*) "concentrations: firstcall, nbq=",nbq
240         write(*,*) "  niq(1:nbq)=",niq(1:nbq)
241         write(*,*) "  aki(1:nbq)=",aki(1:nbq)
242         write(*,*) "  cpi(1:nbq)=",cpi(1:nbq)
243
244         firstcall = .false.
245
246      end if ! if (firstcall)
247
248!     update temperature
249
250      do l = 1,nlayer
251         do ig = 1,ngrid
252            zt(ig,l) = pt(ig,l) + pdt(ig,l)*ptimestep
253         end do
254      end do
255
256!     update tracers
257
258      do l = 1,nlayer
259         do ig = 1,ngrid
260            do i = 1,nbq
261               iq = niq(i)
262               zq(ig,l,iq) = max(1.e-30, pq(ig,l,iq)
263     $                                 + pdq(ig,l,iq)*ptimestep)
264            end do
265         end do
266      end do
267
268!     mmean : mean molecular mass
269!     rnew  : specific gas constant
270
271      mmean(:,:)  = 0.
272
273      do l = 1,nlayer
274         do ig = 1,ngrid
275            do i = 1,nbq
276               iq = niq(i)
277               mmean(ig,l) = mmean(ig,l) + zq(ig,l,iq)/mmol(iq)
278            end do
279            mmean(ig,l) = 1./mmean(ig,l)
280            rnew(ig,l) = 8.314/mmean(ig,l)*1.e3     ! J/kg/K           
281         end do
282      end do
283
284!     cpnew  : specicic heat
285!     akknew : thermal conductivity cofficient
286     
287      cpnew(:,:)  = 0.
288      akknew(:,:) = 0.
289
290      do l = 1,nlayer
291         do ig = 1,ngrid
292            ntot = pplay(ig,l)/(1.381e-23*zt(ig,l))*1.e-6  ! in #/cm3
293            do i = 1,nbq
294               iq = niq(i)
295               ni(iq) = ntot*zq(ig,l,iq)*mmean(ig,l)/mmol(iq)
296               cpnew(ig,l) = cpnew(ig,l) + ni(iq)*cpi(i)
297               akknew(ig,l) = akknew(ig,l) + ni(iq)*aki(i)
298            end do
299            cpnew(ig,l) = cpnew(ig,l)/ntot
300            akknew(ig,l) = akknew(ig,l)/ntot
301         end do
302!        print*, l, mmean(1,l), cpnew(1,l), rnew(1,l)
303      end do
304
305      return
306      end
Note: See TracBrowser for help on using the repository browser.