source: trunk/LMDZ.TITAN/libf/phytitan/calchim.F90 @ 2245

Last change on this file since 2245 was 2099, checked in by jvatant, 6 years ago

Fix a problem of interoperability C-Fortran for picky compilers.
Using iso_c_binding could be a smart future improvement to bring.
--JVO

File size: 22.3 KB
Line 
1SUBROUTINE calchim(ngrid,qy_c,declin,dtchim,            &
2     ctemp,cpphi,cpphis,cplay,cplev,czlay,czlev,dqyc)
3
4  !---------------------------------------------------------------------------------------------------------
5  ! 
6  ! Purpose : Interface subroutine to photochemical C model for Titan GCM.
7  ! -------
8  !           The subroutine computes the chemical processes for a single vertical column.
9  !
10  ! - Only tendencies are returned.
11  ! - With moyzon_ch=.true. and input vectors zonally averaged
12  !   the calculation is done only once per lat. band
13  !
14  ! Authors: + S. Lebonnois : 01/2000 | 09/2003
15  ! -------                  adaptation for Titan 3D : 02/2009
16  !                          adaptation for // : 04/2013
17  !                          extension chemistry up to 1300km : 10/2013
18  !     
19  !          + J. Vatant d'Ollone
20  !               + 02/17 - adaptation for the new generic-forked physics
21  !               + 01/18 - 03/18 - Major transformations :
22  !                   - Upper chemistry fields are now stored in startfi
23  !                     and defined on a pressure grid from Vervack profile
24  !                   - These modifs enables to run chemistry with others resolution than 32x48x55 !
25  !                   - Only the actinic fluxes are still read in a 49-lat input but interp. on lat grid
26  !                   - Chemistry can still be done in 2D
27  !                      -> Calcul. once per band lat and put same tendency in all longi.
28  !                         Check for negs in physiq_mod.
29  !                      -> If procs sharing a lat band, no problem, the calcul will just be done twice.
30  !                      -> Will not work with Dynamico, where the chemistry will have to be done in 3D.
31  !                          ( and there'll be work to do to get rid of averaged fields )
32  !
33  !               + 02/19 : To always have correct photodissociations rates, altitudes sent here by physiq are always
34  !                         calculated with effective g - and with reference to the body not the local surface -
35  !                          even if in physiq we keep altitudes coherent with dynamics !
36  !
37  ! + STILL TO DO : + Replug the interaction with haze (cf titan.old) -> to see with JB.
38  !                 + Use iso_c_binding for the fortran-C exchanges.
39  !---------------------------------------------------------------------------------------------------------
40
41  ! --------------------------------------------------------------------
42  ! Structure :
43  ! -----------
44  !   0.  Declarations
45  !   I.  Init and firstcall
46  !         1. Read and store Vervack profile
47  !         2. Compute planetar averaged atm. properties
48  !         3. Init compound caracteristics
49  !         4. Init photodissociations rates from actinic fluxes
50  !         5. Init chemical reactions
51  !         6. Init eddy diffusion coeff
52  !   II. Loop on latitudes/grid-points
53  !         0. Check on 2D chemistry
54  !           1. Compute atm. properties at grid points
55  !           2. Interpolate photodissociation rates at lat,alt,dec
56  !           3. Read composition
57  !           4. Call main solver gptitan C routine
58  !           5. Calculate output tendencies on advected tracers
59  !           6. Update upper chemistry fields
60  !         0bis. If 2D chemsitry, don't recalculate if needed
61  ! -----------------------------------------------------------------
62
63  USE, INTRINSIC :: iso_c_binding
64  USE comchem_h
65  USE dimphy
66  USE datafile_mod, ONLY: datadir
67  USE comcstfi_mod, ONLY: g, rad, pi, r, kbol
68  USE geometry_mod, ONLY: latitude
69  USE logic_mod, ONLY: moyzon_ch
70  USE moyzon_mod, ONLY: tmoy, playmoy
71
72  IMPLICIT NONE
73
74! ------------------------------------------
75! *********** 0. Declarations *************
76! ------------------------------------------
77
78  ! Arguments
79  ! ---------
80
81  INTEGER, INTENT(IN)                              :: ngrid       ! Number of atmospheric columns.
82  REAL*8, DIMENSION(ngrid,klev,nkim), INTENT(IN)   :: qy_c        ! Chemical species on GCM layers after adv.+diss. (mol/mol).
83  REAL*8, INTENT(IN)                               :: declin      ! Solar declination (rad).
84  REAL*8, INTENT(IN)                               :: dtchim      ! Chemistry timsetep (s).
85  REAL*8, DIMENSION(ngrid,klev),      INTENT(IN)   :: ctemp       ! Mid-layer temperature (K).
86  REAL*8, DIMENSION(ngrid,klev),      INTENT(IN)   :: cpphi       ! Mid-layer geopotential (m2.s-2).
87  REAL*8, DIMENSION(ngrid),           INTENT(IN)   :: cpphis      ! Surface geopotential (m2.s-2).
88  REAL*8, DIMENSION(ngrid,klev),      INTENT(IN)   :: cplay       ! Mid-layer pressure (Pa).
89  REAL*8, DIMENSION(ngrid,klev+1),    INTENT(IN)   :: cplev       ! Inter-layer pressure (Pa).
90  REAL*8, DIMENSION(ngrid,klev),      INTENT(IN)   :: czlay       ! Mid-layer effective altitude (m) : ref = geoid.
91  REAL*8, DIMENSION(ngrid,klev+1),    INTENT(IN)   :: czlev       ! Inter-layer effective altitude (m) ref = geoid.
92
93  REAL*8, DIMENSION(ngrid,klev,nkim), INTENT(OUT)  :: dqyc        ! Chemical species tendencies on GCM layers (mol/mol/s).
94
95  ! Local variables :
96  ! -----------------
97
98  INTEGER :: i , l, ic, ig, igm1
99
100  INTEGER :: dec, idec, ipres, ialt, klat
101
102  REAL*8  :: declin_c  ! Declination (deg).
103  REAL*8  :: factp, factalt, factdec, factlat, krpddec, krpddecp1, krpddecm1
104  REAL*8  :: temp1, logp
105
106  ! Variables sent into chemistry module (must be in double precision)
107  ! ------------------------------------------------------------------
108
109  REAL*8, DIMENSION(nlaykim_tot) :: temp_c  ! Temperature (K).
110  REAL*8, DIMENSION(nlaykim_tot) :: press_c ! Pressure (Pa).
111  REAL*8, DIMENSION(nlaykim_tot) :: phi_c   ! Geopotential (m2.s-2) - actually not sent in chem. module but used to compute alts.
112  REAL*8, DIMENSION(nlaykim_tot) :: nb      ! Density (cm-3).
113
114  REAL*8, DIMENSION(nlaykim_tot,nkim) :: cqy   ! Chemical species in whole column (mol/mol) sent to chem. module.
115  REAL*8, DIMENSION(nlaykim_tot,nkim) :: cqy0  !     "      "     "   "      "        "     before modifs.
116
117  REAL*8  :: surfhaze(nlaykim_tot)
118  REAL*8  :: cprodaer(nlaykim_tot,4), cmaer(nlaykim_tot,4)
119  REAL*8  :: ccsn(nlaykim_tot,4), ccsh(nlaykim_tot,4)
120
121  REAL*8, DIMENSION(nlaykim_tot) :: rmil   ! Mid-layer distance (km) to planetographic center.
122  REAL*8, DIMENSION(nlaykim_tot) :: rinter ! Inter-layer distance (km) to planetographic center (RA grid in chem. module).
123  ! NB : rinter is on nlaykim_tot too, we don't care of the uppermost layer upper boundary altitude.
124
125  ! Saved variables initialized at firstcall
126  ! ----------------------------------------
127
128  LOGICAL, SAVE :: firstcall = .TRUE.
129!$OMP THREADPRIVATE(firstcall)
130
131  REAL*8, DIMENSION(:), ALLOCATABLE, SAVE :: kedd ! Eddy mixing coefficient for upper chemistry (cm^2.s-1)
132!$OMP THREADPRIVATE(kedd)
133
134  REAL*8, DIMENSION(:,:), ALLOCATABLE, SAVE :: md   ! Mean molecular diffusion coefficients (cm^2.s-1)
135  REAL*8, DIMENSION(:),   ALLOCATABLE, SAVE :: mass ! Molar mass of the compounds (g.mol-1)
136!$OMP THREADPRIVATE(mass,md)
137
138  REAL*8, DIMENSION(:), ALLOCATABLE, SAVE :: r1d, ct1d, p1d, t1d ! Vervack profile
139  ! JVO 18 : No threadprivate for those as they'll be read in tcp.ver by master
140
141  REAL*8, DIMENSION(:,:,:,:), ALLOCATABLE, SAVE :: krpd  ! Photodissociations rate table
142  REAL*8, DIMENSION(:,:)    , ALLOCATABLE, SAVE :: krate ! Reactions rate ( photo + chem )
143!$OMP THREADPRIVATE(krpd,krate)
144
145  INTEGER, DIMENSION(:),     ALLOCATABLE, SAVE :: nom_prod, nom_perte
146  INTEGER, DIMENSION(:,:),   ALLOCATABLE, SAVE :: reactif
147  INTEGER, DIMENSION(:,:),   ALLOCATABLE, SAVE :: prod
148  INTEGER, DIMENSION(:,:,:), ALLOCATABLE, SAVE :: perte
149!$OMP THREADPRIVATE(nom_prod,nom_perte,reactif,prod,perte)
150
151  ! TEMPORARY : Dummy parameters without microphysics
152  ! Here to keep the whole stuff running without modif chem. module
153  ! ---------------------------------------------------------------
154
155  INTEGER  :: utilaer(16)
156  INTEGER  :: aerprod = 0
157  INTEGER  :: htoh2   = 0
158
159  ! -----------------------------------------------------------------------
160  ! ***************** I. Initialisations and Firstcall ********************
161  ! -----------------------------------------------------------------------
162
163  IF (firstcall) THEN
164
165     PRINT*, 'CHIMIE, premier appel'
166
167     call check(nlaykim_tot,klev-15,nlrt_kim,nkim)
168
169     ALLOCATE(r1d(131))
170     ALLOCATE(ct1d(131))
171     ALLOCATE(p1d(131))
172     ALLOCATE(t1d(131))
173
174     ALLOCATE(md(nlaykim_tot,nkim))
175     ALLOCATE(mass(nkim))
176     
177     ALLOCATE(kedd(nlaykim_tot))
178     
179     ALLOCATE(krate(nlaykim_tot,nr_kim))
180     ALLOCATE(krpd(nd_kim+1,nlrt_kim,15,nlat_actfluxes))
181
182     ALLOCATE(nom_prod(nkim))
183     ALLOCATE(nom_perte(nkim))
184     ALLOCATE(reactif(5,nr_kim))
185     ALLOCATE(prod(200,nkim))
186     ALLOCATE(perte(2,200,nkim))
187     
188     ! 0. Deal with characters for C-interoperability
189     ! ----------------------------------------------
190     ! NB ( JVO 19 ) : Using iso_c_binding would do things in an even cleaner way !
191     DO ic=1,nkim
192       nomqy_c(ic) = trim(cnames(ic))//char(0) ! Add the C null terminator
193     ENDDO
194     nomqy_c(nkim+1)="HV"//char(0) ! For photodissociations
195
196     ! 1. Read Vervack profile "tcp.ver", once for all
197     ! -----------------------------------------------
198
199!$OMP MASTER
200     OPEN(11,FILE=TRIM(datadir)//'/tcp.ver',STATUS='old')
201     READ(11,*)
202     DO i=1,131
203        READ(11,*) r1d(i), t1d(i), ct1d(i), p1d(i)
204        ! For debug :
205        ! -----------
206        ! PRINT*, "tcp.ver", r1d(i), t1d(i), ct1d(i), p1d(i)
207     ENDDO
208     CLOSE(11)
209!$OMP END MASTER
210!$OMP BARRIER
211
212     ! 2. Calculation of temp_c, densities and altitudes in planetary average
213     ! ----------------------------------------------------------------------
214     
215     ! JVO18 : altitudes are no more calculated in firstcall, as I set kedd in pressure grid
216
217     ! a. For GCM layers we just copy-paste ( assuming that physiq always send correct altitudes ! )
218
219     PRINT*,'Init chemistry : pressure, density, temperature ... :'
220     PRINT*,'level, press_c (mbar), nb (cm-3), temp_c (K)'
221     
222     IF (ngrid.NE.1) THEN
223       DO l=1,klev
224          temp_c(l)  = tmoy(l)                              ! K
225          press_c(l) = playmoy(l)/100.                      ! mbar
226          nb(l)      = 1.e-4*press_c(l) / (kbol*temp_c(l))  ! cm-3
227          PRINT*, l, press_c(l), nb(l), temp_c(l)
228       ENDDO
229     ELSE
230       DO l=1,klev
231          temp_c(l)  = ctemp(1,l)                             ! K
232          press_c(l) = cplay(1,l)/100.                        ! mbar
233          nb(l)      = 1.e-4*press_c(l) / (kbol*temp_c(l))  ! cm-3
234          PRINT*, l, press_c(l), nb(l), temp_c(l)
235       ENDDO
236     ENDIF
237
238     ! b. Extension in upper atmosphere with Vervack profile
239     ! NB : Maybe the transition klev/klev+1 is harsh if T profile different from Vervack ...     
240
241     ipres=1
242     DO l=klev+1,nlaykim_tot
243        press_c(l) = preskim(l-klev) / 100.0
244        DO i=ipres,130
245           IF ( (press_c(l).LE.p1d(i)) .AND. (press_c(l).GT.p1d(i+1)) ) THEN
246              ipres=i
247           ENDIF
248        ENDDO
249        factp = (press_c(l)-p1d(ipres)) / (p1d(ipres+1)-p1d(ipres))
250
251        nb(l)      = exp( log(ct1d(ipres))*(1-factp) + log(ct1d(ipres+1))* factp )
252        temp_c(l)  = t1d(ipres)*(1-factp) + t1d(ipres+1)*factp
253        PRINT*, l , press_c(l), nb(l), temp_c(l)
254     ENDDO
255
256     ! 3. Compounds caracteristics
257     ! ---------------------------
258     mass(:) = 0.0
259     call comp(nomqy_c,nb,temp_c,mass,md)
260     PRINT*,'           Mass'
261     DO ic=1,nkim
262        PRINT*, nomqy_c(ic), mass(ic)
263     ENDDO
264
265     ! 4. Photodissociation rates
266     ! --------------------------
267     call disso(krpd,nlat_actfluxes)
268
269     ! 5. Init. chemical reactions with planetary average T prof.
270     ! ----------------------------------------------------------
271
272     !  NB : Chemical reactions rate are assumed to be constant within the T range of Titan's atm
273     !  so we fill their krate once for all but krate for photodiss will be filled at each timestep
274     
275     call chimie(nomqy_c,nb,temp_c,krate,reactif, &
276          nom_perte,nom_prod,perte,prod)
277
278     ! 6. Eddy mixing coefficients (constant with time and space)
279     ! ----------------------------------------------------------
280     
281     kedd(:) = 1.e3 ! Default value =/= zero
282
283     ! NB : Eddy coeffs (e.g. Lavvas et al 08, Yelle et al 08) in altitude but they're rather linked to pressure
284     !      Below GCM top we have dynamic mixing and for levs < nld=klev-15 the chem. solver ignores diffusion
285
286     !! First calculate kedd for upper chemistry layers
287     !DO l=klev-4,nlaykim_tot
288     !   logp=-log10(press_c(l))
289     !! 2E6 at 400 km ~ 10-2 mbar
290     !   IF     ( logp.ge.2.0 .and. logp.le.3.0 ) THEN
291     !         kedd(l) = 2.e6 * 5.0**(logp-2.0)
292     !! 1E7 at 500 km ~ 10-3 mbar
293     !   ELSE IF     ( logp.ge.3.0 .and. logp.le.4.0 ) THEN
294     !         kedd(l) = 1.e7 * 3.0**(logp-3.0)
295     !! 3E7 above 700 km ~ 10-4 mbar
296     !   ELSEIF ( logp.gt.4.0                   ) THEN
297     !        kedd(l) = 3.e7
298     !   ENDIF
299     !ENDDO
300
301     ! Kedd from (E7) in Vuitton 2019
302     if (ngrid .eq. 1) then ! if 1D no dynamic mixing, we set the kedd in all column
303       DO l=1,nlaykim_tot
304         kedd(l) = 300.0 * ( 1.0E2 / press_c(l) )**1.5 * 3.0E7 /  &
305                 ( 300.0 * ( 1.0E2 / press_c(l) )**1.5 + 3.0E7 )
306       ENDDO
307     else
308       DO l=klev-4,nlaykim_tot
309         ! JVO 18 : We keep the nominal profile in the GCM 5 upper layers
310         !          to have  a correct vertical mixing in the sponge layer
311         kedd(l) = 300.0 * ( 1.0E2 / press_c(l) )**1.5 * 3.0E7 /  &
312               ( 300.0 * ( 1.0E2 / press_c(l) )**1.5 + 3.0E7 )
313       ENDDO
314     endif
315     
316     if (ngrid .gt. 1) then ! not in 1D, no dynamic mixing
317       ! Then adjust 10 layers profile fading to default value depending on kedd(ptop)
318       DO l=klev-15,klev-5
319          temp1   = ( log10(press_c(l)/press_c(klev-15)) ) / ( log10(press_c(klev-4)/press_c(klev-15)) )
320          kedd(l) = 10.**( 3.0 + log10(kedd(klev-4)/1.e3) * temp1 )
321       ENDDO
322     endif
323 
324     firstcall = .FALSE.
325  ENDIF  ! firstcall
326 
327  declin_c = declin*180./pi
328
329  ! -----------------------------------------------------------------------
330  ! *********************** II. Loop on latitudes *************************
331  ! -----------------------------------------------------------------------
332 
333  DO ig=1,ngrid
334
335    IF (ig.eq.1) THEN
336        igm1=1
337     ELSE
338        igm1=ig-1
339     ENDIF
340
341     ! If 2D chemistry, trick to do the calculation only once per latitude band within the chunk
342     ! NB1 : Will be obsolete with DYNAMICO, the chemistry will necessarly be 3D
343     ! NB2 : Test of same latitude with dlat=0.1 : I think that if you run sims better than 1/10th degree then
344     ! either it's with Dynamico and doesn't apply OR it is more than enough in terms of "preco / calc time" !
345     ! -------------------------------------------------------------------------------------------------------
346
347     IF ( ( moyzon_ch .AND. ( ig.EQ.1 .OR. (ABS(latitude(ig)-latitude(igm1)).GT.0.1*pi/180.0)) ) .OR. (.NOT. moyzon_ch) ) THEN
348
349        ! 1. Compute altitude for the grid point with hydrostat. equilib.
350        ! ---------------------------------------------------------------
351
352        ! a. For GCM layers we just copy-paste
353        ! JVO 19 : Now physiq always sent correct altitudes with effective g for chemistry ( even if it's not the case in physiq )
354
355        DO l=1,klev
356           rinter(l)  = (czlev(ig,l)+rad)/1000.0             ! km
357           rmil(l)    = (czlay(ig,l)+rad)/1000.0             ! km
358           temp_c(l)  = ctemp(ig,l)                          ! K
359           phi_c(l)   = cpphi(ig,l)                          ! m2.s-2
360           press_c(l) = cplay(ig,l)/100.                     ! mbar
361           nb(l)      = 1.e-4*press_c(l) / (kbol*temp_c(l))  ! cm-3
362        ENDDO
363        rinter(klev+1)=(czlev(ig,klev+1)+rad)/1000.
364
365        ! b. Extension in upper atmosphere with Vervack profile
366
367        ipres=1
368        DO l=klev+1,nlaykim_tot
369          press_c(l) = preskim(l-klev) / 100.0
370          DO i=ipres,130
371              IF ( (press_c(l).LE.p1d(i)) .AND. (press_c(l).GT.p1d(i+1)) ) THEN
372                ipres=i
373              ENDIF
374          ENDDO
375          factp = (press_c(l)-p1d(ipres)) / (p1d(ipres+1)-p1d(ipres))
376
377          nb(l)      = exp( log(ct1d(ipres))*(1-factp) + log(ct1d(ipres+1))* factp )
378          temp_c(l)  = t1d(ipres)*(1-factp) + t1d(ipres+1)*factp
379        ENDDO
380 
381        ! We build altitude with hydrostatic equilibrium on preskim grid with Vervack profile
382        ! ( keeping in mind that preskim is built based on Vervack profile with dz=10km )
383
384        DO l=klev+1,nlaykim_tot
385
386           ! Compute geopotential on the upper grid with effective g to have correct altitudes
387           ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
388
389           temp1 = 0.5*(temp_c(l-1)+temp_c(l)) ! interlayer temp
390           phi_c(l) = phi_c(l-1) + r*temp1*log(press_c(l-1)/press_c(l)) ! Geopotential assuming hydrostatic equilibrium
391
392           rmil(l)  =  ( g*rad*rad / (g*rad - ( phi_c(l) + cpphis(ig) ) ) ) / 1000.0 ! z(phi) with g varying with altitude with reference to the geoid
393        ENDDO
394
395        DO l=klev+2,nlaykim_tot
396          rinter(l) = 0.5*(rmil(l-1) + rmil(l)) ! should be balanced with the intermediate pressure rather than 0.5
397        ENDDO
398
399        ! 2. From krpd, compute krate for dissociations (declination-latitude-altitude interpolation)
400        ! -------------------------------------------------------------------------------------------
401
402        ! a. Calculate declination dependence
403
404        if ((declin_c*10+267).lt.14.) then
405           idec = 0
406           dec  = 0
407        else
408           if ((declin_c*10+267).gt.520.) then
409              idec = 14
410              dec  = 534
411           else
412              idec = 1
413              dec  = 27
414              do while( (declin_c*10+267).ge.real(dec+20) )
415                 dec  = dec+40
416                 idec = idec+1
417              enddo
418           endif
419        endif
420        if ((declin_c.ge.-24.).and.(declin_c.le.24.)) then
421           factdec = ( declin_c - (dec-267)/10. ) / 4.
422        else
423           factdec = ( declin_c - (dec-267)/10. ) / 2.7
424        endif
425
426        ! b. Calculate klat for interpolation on fixed latitudes of actinic fluxes input
427
428        klat=1
429        DO i=1,nlat_actfluxes
430          IF (latitude(ig).LT.lat_actfluxes(i)) klat=i
431        ENDDO
432        IF (klat==nlat_actfluxes) THEN ! avoid rounding problems
433          klat    = nlat_actfluxes-1
434          factlat = 1.0
435        ELSE
436          factlat = (latitude(ig)-lat_actfluxes(klat))/(lat_actfluxes(klat+1)-lat_actfluxes(klat))
437        ENDIF
438
439        ! c. Altitude loop
440
441        DO l=1,nlaykim_tot
442
443           ! Calculate ialt for interpolation in altitude (krpd every 2 km)
444           ialt    = int((rmil(l)-rad/1000.)/2.)+1
445           factalt = (rmil(l)-rad/1000.)/2.-(ialt-1)
446
447           ! Altitude can go above top limit of UV levels - in this case we keep the 1310km top fluxes
448           IF (ialt.GT.nlrt_kim-1) THEN
449             ialt     = nlrt_kim-1 ! avoid out-of-bound array
450             factalt  = 1.0
451           ENDIF
452
453           DO i=1,nd_kim+1 ! nd_kim+1 is dissociation of N2 by GCR
454
455                 krpddec   =   (   krpd(i,ialt  ,idec+1,klat)   * (1.0-factalt)                   &
456                                 + krpd(i,ialt+1,idec+1,klat)   * factalt       ) * (1.0-factlat) &
457                             + (   krpd(i,ialt  ,idec+1,klat+1) * (1.0-factalt)                   &
458                                 + krpd(i,ialt+1,idec+1,klat+1) * factalt       ) * factlat
459
460              if      ( factdec.lt.0. ) then
461                 krpddecm1 =   (   krpd(i,ialt  ,idec  ,klat)   * (1.0-factalt)                   &
462                                 + krpd(i,ialt+1,idec  ,klat)   * factalt       ) * (1.0-factlat) &
463                             + (   krpd(i,ialt  ,idec  ,klat+1) * (1.0-factalt)                   &
464                                 + krpd(i,ialt+1,idec  ,klat+1) * factalt       ) * factlat
465                 krate(l,i) = krpddecm1 * abs(factdec) + krpddec   * ( 1.0 + factdec)
466              else if ( factdec.gt.0. ) then
467                 krpddecp1 =   (   krpd(i,ialt  ,idec+2,klat)   * (1.0-factalt)                   &
468                                 + krpd(i,ialt+1,idec+2,klat)   * factalt       ) * (1.0-factlat) &
469                             + (   krpd(i,ialt  ,idec+2,klat+1) * (1.0-factalt)                   &
470                                 + krpd(i,ialt+1,idec+2,klat+1) * factalt       ) * factlat
471                 krate(l,i) = krpddecp1 * factdec      + krpddec   * ( 1.0 - factdec)
472              else if ( factdec.eq.0. ) then
473                 krate(l,i) = krpddec
474              endif
475
476           ENDDO ! i=1,nd_kim+1
477        ENDDO ! l=1,nlaykim_tot
478
479        ! 3. Read composition
480        ! -------------------
481
482        DO ic=1,nkim
483           DO l=1,klev
484              cqy(l,ic) = qy_c(ig,l,ic) ! advected tracers for the GCM part converted to molar frac.
485           ENDDO
486           
487           DO l=1,nlaykim_up
488              cqy(klev+l,ic) = ykim_up(ic,ig,l) ! ykim_up for the upper atm.
489           ENDDO
490        ENDDO
491
492        cqy0(:,:) = cqy(:,:) ! Stores compo. before modifs
493
494        ! 4. Call main Titan chemistry C routine
495        ! --------------------------------------
496
497        call gptitan(rinter,temp_c,nb,                  &
498             nomqy_c,cqy,                               &
499             dtchim,latitude(ig)*180./pi,mass,md,       &
500             kedd,krate,reactif,                        &
501             nom_prod,nom_perte,prod,perte,             &
502             aerprod,utilaer,cmaer,cprodaer,ccsn,ccsh,  &
503             htoh2,surfhaze)
504
505        ! 5. Calculates tendencies on composition for advected tracers
506        ! ------------------------------------------------------------
507        DO ic=1,nkim
508           DO l=1,klev
509              dqyc(ig,l,ic) = (cqy(l,ic) - cqy0(l,ic))/dtchim ! (mol/mol/s)
510           ENDDO
511        ENDDO
512
513        ! 6. Update ykim_up
514        ! -----------------
515        DO ic=1,nkim
516           DO l=1,nlaykim_up
517              ykim_up(ic,ig,l) = cqy(klev+l,ic)
518           ENDDO
519        ENDDO
520        ! NB: The full vertical composition grid will be created only for the outputs
521
522
523     ELSE ! In 2D chemistry, if following grid point at same latitude, same zonal mean so don't do calculations again !
524        dqyc(ig,:,:)    = dqyc(igm1,:,:) ! will be put back in 3D with longitudinal variations assuming same relative tendencies within a lat band
525        ykim_up(:,ig,:) = ykim_up(:,igm1,:) ! no horizontal mixing in upper layers -> no longitudinal variations
526     ENDIF
527
528  ENDDO
529
530END SUBROUTINE calchim
Note: See TracBrowser for help on using the repository browser.