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

Last change on this file since 3026 was 2368, checked in by jvatant, 5 years ago

Titan GCM:Fix missing MESOSCALE cpp key in calchim

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