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

Last change on this file since 1680 was 1672, checked in by jvatant, 8 years ago

Re-plug chemistry of old Titan
+ minor correction on reading Titan's startfiles and haze routine
JVO

File size: 16.4 KB
Line 
1      SUBROUTINE calchim(nlon,ny,qy_c,nomqy_c,declin_rad,ls_rad, &
2                         dtchim,ctemp,cplay,cplev,czlay,czlev,dqyc,NLEV)
3     
4!--------------------------------------------------------------------
5
6!     Introduction d une routine chimique
7!
8!     Auteurs: S. Lebonnois,  01/2000 | 09/2003
9!            adaptation pour Titan 3D: 02/2009
10!            adaptation pour // : 04/2013
11!            extension chimie jusqu a 1300km : 10/2013
12!     
13!              J. Vatant d'Ollone, 02/2017
14!               + adaptation pour le nouveau titan issu du generic
15!
16!---------------------------------------------------------------------
17!
18     
19      use dimphy
20     
21      use datafile_mod, only: datadir
22     
23      use comcstfi_mod, only: rad, pi, kbol
24      use moyzon_mod, only: tmoy,playmoy,zlaymoy,zlevmoy,klat
25      use mod_grid_phy_lmdz, only: nbp_lat
26     
27      implicit none
28
29!   Parameters ( dans common_mod in old titan)
30!              + doivent etre en accord avec titan.h
31!   -------------------------------------------------
32
33      INTEGER,PARAMETER ::   NC   = 44     ! Nb de composes dans la chimie
34      INTEGER,PARAMETER ::   ND   = 54     ! Nb de photodissociations
35      INTEGER,PARAMETER ::   NR   = 377    ! Nb de reactions dans la chimie
36      INTEGER,PARAMETER ::   NLRT = 650    ! Pour l'UV 650 niveaux de 2 km
37     
38!     Dummy parameters without microphysics
39!      + Just here to keep the whole stuff running without modif C sources
40!     ---------------------------------------------------------------------
41
42      INTEGER  :: utilaer(16)
43      INTEGER  :: aerprod = 0
44      INTEGER  :: htoh2   = 0
45
46!    Arguments
47!    ---------
48
49      INTEGER      nlon                   ! nb of horiz points
50      INTEGER      ny                     ! nb de composes (nqmax-nmicro)
51      REAL*8         qy_c(nlon,klev,NC)     ! Especes chimiques apres adv.+diss.
52      character*10 nomqy_c(NC+1)          ! Noms des especes chimiques
53      REAL*8         declin_rad,ls_rad      ! declinaison et long solaire en radians
54      REAL*8         dtchim                 ! pas de temps chimie
55      REAL*8         ctemp(nlon,klev)       ! Temperature
56      REAL*8         cplay(nlon,klev)       ! pression (Pa)
57      REAL*8         cplev(nlon,klev+1)     ! pression intercouches (Pa)
58      REAL*8         czlay(nlon,klev)       ! altitude (m)
59      REAL*8         czlev(nlon,klev+1)     ! altitude intercouches (m)
60     
61      REAL*8         dqyc(nlon,klev,NC)     ! Tendances especes chimiques
62     
63      INTEGER      NLEV                   ! nbp_lev+70 - a changer si =/ 55 layers ?
64     
65!    Local variables :
66!    -----------------
67
68      integer i,j,l,ic,jm1
69
70! variables envoyees dans la chimie: double precision
71
72      REAL*8  temp_c(NLEV)
73      REAL*8  press_c(NLEV)   ! T,p(mbar) a 1 lat donnee
74      REAL*8  cqy(NLEV,NC)    ! frac mol qui seront modifiees
75      REAL*8  cqy0(NLEV,NC)    ! frac mol avant modif
76      REAL*8  surfhaze(NLEV)
77      REAL*8  cprodaer(NLEV,4),cmaer(NLEV,4)
78      REAL*8  ccsn(NLEV,4),ccsh(NLEV,4)
79! rmil: milieu de couche, grille pour K,D,p,ct,T,y
80! rinter: intercouche (grille RA dans la chimie)
81      REAL*8  rmil(NLEV),rinter(NLEV),nb(NLEV)
82      REAL*8,save,allocatable :: kedd(:)
83
84      REAL*8  a,b
85      character str1*1,str2*2
86      REAL*8  latit
87      character*20 formt1,formt2,formt3
88     
89!     variables locales initialisees au premier appel
90
91      LOGICAL firstcall
92      DATA    firstcall/.true./
93      SAVE    firstcall
94     
95      integer dec,declinint,ialt
96      REAL*8  declin_c                       ! declinaison en degres
97      real*8  factalt,factdec,krpddec,krpddecp1,krpddecm1
98      real*8  duree
99      REAL*8,save :: mass(NC)
100      REAL*8,save,allocatable :: md(:,:)
101      REAL*8,save :: botCH4
102      DATA  botCH4/0.05/
103      real*8,save ::  r1d(131),ct1d(131),p1d(131),t1d(131) ! lecture tcp.ver
104      REAL*8,save,allocatable :: krpd(:,:,:,:),krate(:,:)
105      integer,save :: reactif(5,NR),nom_prod(NC),nom_perte(NC)
106      integer,save :: prod(200,NC),perte(2,200,NC)
107      character  dummy*30,fich*7,ficdec*15,curdec*15,name*10
108      real*8  ficalt,oldq,newq,zalt
109      logical okfic
110
111
112!-----------------------------------------------------------------------
113!***********************************************************************
114!
115!    Initialisations :
116!    ----------------
117
118      duree = dtchim  ! passage en real*4 pour appel a routines C
119
120      IF (firstcall) THEN
121           
122        print*,'CHIMIE, premier appel'
123       
124! ************************************
125! Au premier appel, initialisation de toutes les variables
126! pour les routines de la chimie.
127! ************************************
128
129        allocate(kedd(NLEV))
130
131        allocate(krpd(15,ND+1,NLRT,nbp_lat),krate(NLEV,NR),md(NLEV,NC))
132
133! Verification du nombre de composes: coherence common_mod et nqmax-nmicro
134! ----------------------------------
135
136        if (ny.ne.NC) then
137           print*,'PROBLEME de coherence dans le nombre de composes:',ny,NC
138           STOP
139        endif
140
141! calcul de temp_c, densites et press_c en moyenne planetaire :
142! --------------------------------------------------------------
143
144        print*,'pression, densites et temp (init chimie):'
145        print*,'level, press_c, nb, temp_c'
146        DO l=1,klev
147         rinter(l)  = (zlevmoy(l)+rad)/1000.
148         rmil(l)    = (zlaymoy(l)+rad)/1000.
149!     temp_c (K):
150         temp_c(l)  = tmoy(l)
151!     press_c (mbar):
152         press_c(l) = playmoy(l)/100.
153!     nb (cm-3):
154         nb(l) = 1.e-4*press_c(l) / (kbol*temp_c(l))
155         print*,l,rmil(l)-rad/1000.,press_c(l),nb(l),temp_c(l)
156        ENDDO
157        rinter(klev+1)=(zlevmoy(klev+1)+rad)/1000.
158
159! au-dessus du GCM, dr regulier et rinter(NLEV)=1290+2575 km.
160       do l=klev+2,NLEV
161         rinter(l) = rinter(klev+1) &
162               + (l-klev-1)*(3865.-rinter(klev+1))/(NLEV-klev-1)
163         rmil(l-1) = (rinter(l-1)+rinter(l))/2.
164       enddo
165       rmil(NLEV) = rinter(NLEV)+(rinter(NLEV)-rinter(NLEV-1))/2.
166
167! lecture de tcp.ver, une seule fois
168! remplissage r1d,t1d,ct1d,p1d
169        open(11,file=TRIM(datadir)//'/tcp.ver',status='old')
170        read(11,*) dummy
171        do i=1,131
172          read(11,*) r1d(i),t1d(i),ct1d(i),p1d(i)
173          !print*,"TCP.VER ",r1d(i),t1d(i),ct1d(i),p1d(i)
174        enddo
175        close(11)
176
177! extension pour klev+1 a NLEV avec tcp.ver
178! la jonction klev/klev+1 est brutale... Tant pis ?
179        ialt=1
180        do l=klev+1,NLEV
181           zalt=rmil(l)-rad/1000.
182           do i=ialt,130
183            if ((zalt.ge.r1d(i)).and.(zalt.lt.r1d(i+1))) then
184              ialt=i
185            endif
186           enddo
187           factalt = (zalt-r1d(ialt))/(r1d(ialt+1)-r1d(ialt))
188           press_c(l) = exp(  log(p1d(ialt))   * (1-factalt)  &
189                           + log(p1d(ialt+1)) * factalt    )
190           nb(l)      = exp(  log(ct1d(ialt))   * (1-factalt) & 
191                           + log(ct1d(ialt+1)) * factalt    )
192           temp_c(l)  = t1d(ialt) * (1-factalt) + t1d(ialt+1) * factalt
193!          print*,l,zalt,press_c(l),nb(l),temp_c(l)
194        enddo
195       
196! caracteristiques des composes:       
197! -----------------------------
198        mass(:) = 0.0
199        call comp(nomqy_c,nb,temp_c,mass,md)
200        print*,'           Mass'
201        do ic=1,NC
202          print*,nomqy_c(ic),mass(ic)
203!         if(nomqy_c(ic).eq.'CH4') print*,"MD=",md(:,ic)
204        enddo   
205               
206! taux de photodissociations:
207! --------------------------
208        call disso(krpd,nbp_lat)
209
210! reactions chimiques:
211! -------------------
212        call chimie(nomqy_c,nb,temp_c,krate,reactif, &
213                    nom_perte,nom_prod,perte,prod)
214!        print*,'nom_prod, nom_perte:'
215!        do ic=1,NC
216!          print*,nom_prod(ic),nom_perte(ic)
217!        enddo
218!        print*,'premieres prod, perte(1:reaction,2:compagnon):'
219!        do ic=1,NC
220!          print*,prod(1,ic),perte(1,1,ic),perte(2,1,ic)
221!        enddo
222
223!       l = klev-3
224!       print*,'krate a p=',press_c(l),' reactifs et produits:'
225!       do ic=1,ND+1
226!         print*,ic,krpd(7,ic,l,4)*nb(l),"  ",
227!    .     nomqy_c(reactif(1,ic)+1),
228!    .     nomqy_c(reactif(2,ic)+1),nomqy_c(reactif(3,ic)+1),
229!    .     nomqy_c(reactif(4,ic)+1),nomqy_c(reactif(5,ic)+1)
230!       enddo
231!       do ic=ND+2,NR
232!         print*,ic,krate(l,ic),"  ",
233!    .     nomqy_c(reactif(1,ic)+1),
234!    .     nomqy_c(reactif(2,ic)+1),nomqy_c(reactif(3,ic)+1),
235!    .     nomqy_c(reactif(4,ic)+1),nomqy_c(reactif(5,ic)+1)
236!       enddo
237
238
239!   init kedd
240!   ---------
241!   kedd stays constant with time and space
242!   => linked to pressure rather than altitude...
243 
244      kedd(:) = 5e2 ! valeur mise par defaut 
245               ! UNITE ? doit etre ok pour gptitan
246      do l=1,NLEV
247       zalt=rmil(l)-rad/1000.  ! z en km
248       if     ((zalt.ge.250.).and.(zalt.le.400.)) then
249         kedd(l) = 10.**(3.+(zalt-250.)/50.)
250! 1E3 at 250 km
251! 1E6 at 400 km
252       elseif ((zalt.gt.400.).and.(zalt.le.600.)) then
253         kedd(l) = 10.**(6.+1.3*(zalt-400.)/200.)
254! 2E7 at 600 km
255       elseif ((zalt.gt.600.).and.(zalt.le.900.)) then
256         kedd(l) = 10.**(7.3+0.7*(zalt-600.)/300.)
257! 1E8 above 900 km
258       elseif ( zalt.gt.900.                    ) then
259        kedd(l) = 1.e8
260       endif
261      enddo
262
263      ENDIF  ! premier appel
264
265!***********************************************************************
266!-----------------------------------------------------------------------
267
268!   calcul declin_c (en degres)
269!   ---------------------------
270
271      declin_c = declin_rad*180./pi
272!      print*,'declinaison en degre=',declin_c
273       
274!***********************************************************************
275!***********************************************************************
276!
277!                BOUCLE SUR LES LATITUDES
278!
279      DO j=1,nlon
280     
281      if (j.eq.1) then
282         jm1=1
283      else
284         jm1=j-1
285      endif
286
287      if((j.eq.1).or.(klat(j).ne.klat(jm1))) then
288
289!***********************************************************************
290!***********************************************************************
291
292!-----------------------------------------------------------------------
293!
294!   Distance radiale (intercouches, en km)
295!   ----------------------------------------
296
297       do l=1,klev
298         rinter(l) = (rad+czlev(j,l))/1000.
299         rmil(l)   = (rad+czlay(j,l))/1000.
300!        print*,rinter(l)
301       enddo
302       rinter(klev+1)=(rad+czlev(j,klev+1))/1000.
303
304! au-dessus du GCM, dr regulier et rinter(NLEV)=1290+2575 km.
305       do l=klev+2,NLEV
306         rinter(l) = rinter(klev+1)  &
307               + (l-klev-1)*(3865.-rinter(klev+1))/(NLEV-klev-1)
308         rmil(l-1) = (rinter(l-1)+rinter(l))/2.
309       enddo
310       rmil(NLEV) = rinter(NLEV)+(rinter(NLEV)-rinter(NLEV-1))/2.
311
312!-----------------------------------------------------------------------
313!
314!   Temperature, pression (mbar), densite (cm-3)
315!   -------------------------------------------
316
317       DO l=1,klev
318!     temp_c (K):
319               temp_c(l)  = ctemp(j,l)
320!     press_c (mbar):
321               press_c(l) = cplay(j,l)/100.
322!     nb (cm-3):
323               nb(l) = 1.e-4*press_c(l) / (kbol*temp_c(l))
324       ENDDO
325! extension pour klev+1 a NLEV avec tcp.ver
326       ialt=1
327       do l=klev+1,NLEV
328           zalt=rmil(l)-rad/1000.
329           do i=ialt,130
330            if ((zalt.ge.r1d(i)).and.(zalt.lt.r1d(i+1))) then
331              ialt=i
332            endif
333           enddo
334           factalt = (zalt-r1d(ialt))/(r1d(ialt+1)-r1d(ialt))
335           press_c(l) = exp(  log(p1d(ialt))   * (1-factalt) &
336                           + log(p1d(ialt+1)) * factalt    )
337           nb(l)      = exp(  log(ct1d(ialt))   * (1-factalt) &
338                           + log(ct1d(ialt+1)) * factalt    )
339           temp_c(l)  = t1d(ialt) * (1-factalt) + t1d(ialt+1) * factalt
340       enddo
341               
342!-----------------------------------------------------------------------
343!
344!   passage krpd => krate
345!   ---------------------
346! initialisation krate pour dissociations
347
348      if ((declin_c*10+267).lt.14.) then
349          declinint = 0
350          dec       = 0
351      else
352       if ((declin_c*10+267).gt.520.) then
353          declinint = 14
354          dec       = 534
355       else
356          declinint = 1
357          dec       = 27
358          do while( (declin_c*10+267).ge.real(dec+20) )
359            dec       = dec+40
360            declinint = declinint+1
361          enddo
362       endif
363      endif
364      if ((declin_c.ge.-24.).and.(declin_c.le.24.)) then
365          factdec = ( declin_c - (dec-267)/10. ) / 4.
366      else
367          factdec = ( declin_c - (dec-267)/10. ) / 2.7
368      endif
369
370      do l=1,NLEV
371
372! INTERPOL EN ALT POUR k (krpd tous les deux km)
373        ialt    = int((rmil(l)-rad/1000.)/2.)+1
374        factalt = (rmil(l)-rad/1000.)/2.-(ialt-1)
375
376        do i=1,ND+1
377          krpddecm1 = krpd(declinint  ,i,ialt  ,klat(j)) * (1-factalt) &
378                   + krpd(declinint  ,i,ialt+1,klat(j)) * factalt
379          krpddec   = krpd(declinint+1,i,ialt  ,klat(j)) * (1-factalt) &
380                   + krpd(declinint+1,i,ialt+1,klat(j)) * factalt
381          krpddecp1 = krpd(declinint+2,i,ialt  ,klat(j)) * (1-factalt) &
382                   + krpd(declinint+2,i,ialt+1,klat(j)) * factalt
383
384                    ! ND+1 correspond a la dissociation de N2 par les GCR
385          if ( factdec.lt.0. ) then
386        krate(l,i) = krpddecm1 * abs(factdec)  &
387                  + krpddec   * ( 1 + factdec)
388          endif
389          if ( factdec.gt.0. ) then
390        krate(l,i) = krpddecp1 * factdec       &
391                  + krpddec   * ( 1 - factdec)
392          endif
393          if ( factdec.eq.0. ) krate(l,i) = krpddec
394        enddo       
395      enddo       
396
397!-----------------------------------------------------------------------
398!
399!   composition
400!   ------------
401
402       do ic=1,NC
403        do l=1,klev
404          cqy(l,ic)  = qy_c(j,l,ic)
405          cqy0(l,ic) = cqy(l,ic)
406        enddo
407       enddo
408
409! lecture du fichier compo_klat(j) (01 à 49) pour klev+1 a NLEV
410
411      WRITE(str2,'(i2.2)') klat(j)
412      fich = "comp_"//str2//".dat"
413      inquire (file=fich,exist=okfic)
414      if (okfic) then
415       open(11,file=fich,status='old')
416! premiere ligne=declin
417       read(11,'(A15)') ficdec
418       write(curdec,'(E15.5)') declin_c
419! si la declin est la meme: ce fichier a deja ete reecrit
420! on lit la colonne t-1 au lieu de la colonne t
421! (cas d une bande de latitude partagee par 2 procs)
422       do ic=1,NC
423        read(11,'(A10)') name
424        if (name.ne.nomqy_c(ic)) then
425          print*,"probleme lecture ",fich
426          print*,name,nomqy_c(ic)
427          stop
428        endif
429        if (ficdec.eq.curdec) then
430          do l=klev+1,NLEV
431            read(11,*) ficalt,cqy(l,ic),newq
432          enddo
433        else
434          do l=klev+1,NLEV
435            read(11,*) ficalt,oldq,cqy(l,ic)
436          enddo
437        endif
438       enddo
439       close(11)
440      else       ! le fichier n'est pas la
441       do ic=1,NC
442        do l=klev+1,NLEV
443          cqy(l,ic)=cqy(klev,ic)
444        enddo
445       enddo
446      endif
447      cqy0 = cqy
448
449!-----------------------------------------------------------------------
450!
451!   Appel de chimietitan
452!   --------------------
453       
454       call gptitan(rinter,temp_c,nb,                          &
455                   nomqy_c,cqy,                               &
456                   duree,(klat(j)-1),mass,md,                 &
457                   kedd,botCH4,krate,reactif,                 &
458                   nom_prod,nom_perte,prod,perte,             &
459                   aerprod,utilaer,cmaer,cprodaer,ccsn,ccsh,  &
460                   htoh2,surfhaze)
461       
462!   Tendances composition
463!   ---------------------
464     
465       do ic=1,NC
466        do l=1,klev
467          dqyc(j,l,ic) = (cqy(l,ic) - cqy0(l,ic))/dtchim  ! en /s
468        enddo
469       enddo
470
471
472!-----------------------------------------------------------------------
473!
474!   sauvegarde compo sur NLEV
475!   -------------------------
476
477! dans fichier compo_klat(j) (01 à 49)
478       
479      open(11,file=fich,status='replace')
480! premiere ligne=declin
481      write(11,'(E15.5)') declin_c
482      do ic=1,NC
483       write(11,'(A10)') nomqy_c(ic)
484       do l=klev+1,NLEV
485        write(11,'(E15.5,1X,E15.5,1X,E15.5)') rmil(l)-rad/1000.,cqy0(l,ic),cqy(l,ic)
486       enddo
487      enddo
488      close(11)
489
490!***********************************************************************
491!***********************************************************************
492!              FIN: BOUCLE SUR LES LATITUDES
493
494      else      ! same latitude, we don't do calculations again
495        dqyc(j,:,:) = dqyc(jm1,:,:)
496      endif
497
498      ENDDO
499     
500!***********************************************************************
501!***********************************************************************
502
503      firstcall = .false.
504     
505      end subroutine calchim
Note: See TracBrowser for help on using the repository browser.