source: trunk/LMDZ.MARS/libf/aeronomars/moldiff_red.F90 @ 2599

Last change on this file since 2599 was 2467, checked in by emillour, 4 years ago

Mars GCM:

  • bug fix in the diagnostic of H, H2, and D total escape which was wrong when running in parallel. The escape rates can now also be included as outputs in diagfi.

FGG

File size: 40.3 KB
Line 
1subroutine moldiff_red(ngrid,nlayer,nq,pplay,pplev,pt,pdt,pq,pdq,&
2                       ptimestep,zzlay,pdteuv,pdtconduc,pdqdiff,&
3                       PhiEscH,PhiEscH2,PhiEscD)
4
5use tracer_mod, only: noms, mmol
6use geometry_mod, only: cell_area
7use planetwide_mod, only: planetwide_sumval
8
9implicit none
10
11include "diffusion.h"
12
13! July 2014 JYC ADD BALISTIC Transport coupling to compute wup for H and H2
14
15
16
17!
18! Input/Output
19!
20      integer,intent(in) :: ngrid ! number of atmospheric columns
21      integer,intent(in) :: nlayer ! number of atmospheric layers
22      integer,intent(in) :: nq ! number of advected tracers
23      real ptimestep
24      real pplay(ngrid,nlayer)
25      real zzlay(ngrid,nlayer)
26      real pplev(ngrid,nlayer+1)
27      real pq(ngrid,nlayer,nq)
28      real pdq(ngrid,nlayer,nq)
29      real pt(ngrid,nlayer)
30      real pdt(ngrid,nlayer)
31      real pdteuv(ngrid,nlayer)
32      real pdtconduc(ngrid,nlayer)
33      real pdqdiff(ngrid,nlayer,nq)
34      real*8 PhiEscH,PhiEscH2,PhiEscD
35!
36! Local
37!
38     
39
40!      real hco2(ncompdiff),ho
41
42      integer,dimension(nq) :: indic_diff
43      integer ig,iq,nz,l,k,n,nn,p,ij0
44      integer istep,il,gcn,ntime,nlraf
45      real*8 masse
46      real*8 masseU,kBolt,RgazP,Rmars,Grav,Mmars
47      real*8 rho0,D0,T0,H0,time0,dZ,time,dZraf,tdiff,Zmin,Zmax
48      real*8 FacEsc,invsgmu
49      real*8 PhiauxH(ngrid),PhiauxH2(ngrid),PhiauxD(ngrid)
50      real*8 hp(nlayer)
51      real*8 pp(nlayer)
52      real*8 pint(nlayer)
53      real*8 tt(nlayer),tnew(nlayer),tint(nlayer)
54      real*8 zz(nlayer)
55      real*8,dimension(:,:),allocatable,save :: qq,qnew,qint,FacMass
56      real*8,dimension(:,:),allocatable,save :: rhoK,rhokinit
57      real*8 rhoT(nlayer)
58      real*8 dmmeandz(nlayer)
59      real*8 massemoy(nlayer)
60      real*8,dimension(:),allocatable :: Praf,Traf,Rraf,Mraf,Nraf,Draf,Hraf,Wraf
61      real*8,dimension(:),allocatable :: Zraf,Tdiffraf
62      real*8,dimension(:),allocatable :: Prafold,Mrafold
63      real*8,dimension(:,:),allocatable :: Qraf,Rrafk,Nrafk
64      real*8,dimension(:,:),allocatable :: Rrafkold
65      real*8,dimension(:,:),allocatable :: Drafmol,Hrafmol,Wrafmol,Tdiffrafmol
66      real*8,dimension(:),allocatable :: Atri,Btri,Ctri,Dtri,Xtri,Tad,Dad,Zad,rhoad
67      real*8,dimension(:),allocatable :: alpha,beta,gama,delta,eps
68      real*8,dimension(:),allocatable,save ::  wi,Wad,Uthermal,Lambdaexo,Hspecie
69      real*8,dimension(:),allocatable,save :: Mtot1,Mtot2,Mraf1,Mraf2
70      integer,parameter :: ListeDiffNb=16
71      character(len=20),dimension(ListeDiffNb) :: ListeDiff
72      real*8,parameter :: pi=3.141592653
73      real*8,parameter :: g=3.72d0
74!ccccccccccccccccccccccccccccccccccccccccccccccccccccccc
75!     tracer numbering in the molecular diffusion
76!ccccccccccccccccccccccccccccccccccccccccccccccccccccccc
77!  We need the index of escaping species
78
79      integer,save :: i_h2 
80      integer,save :: i_h
81      integer,save :: i_d
82! vertical index limit for the molecular diffusion
83      integer,save :: il0 
84
85! Tracer indexes in the GCM:
86!      integer,save :: g_co2=0
87!      integer,save :: g_co=0
88!      integer,save :: g_o=0
89!      integer,save :: g_o1d=0
90!      integer,save :: g_o2=0
91!      integer,save :: g_o3=0
92!      integer,save :: g_h=0
93!      integer,save :: g_h2=0
94!      integer,save :: g_oh=0
95!      integer,save :: g_ho2=0
96!      integer,save :: g_h2o2=0
97!      integer,save :: g_n2=0
98!      integer,save :: g_ar=0
99!      integer,save :: g_h2o=0
100!      integer,save :: g_n=0
101!      integer,save :: g_no=0
102!      integer,save :: g_no2=0
103!      integer,save :: g_n2d=0
104!      integer,save :: g_oplus=0
105!      integer,save :: g_hplus=0
106
107      integer,save :: ncompdiff
108      integer,dimension(:),allocatable,save :: gcmind ! array of GCM indexes
109      integer ierr,compteur
110
111      logical,save :: firstcall=.true.
112!      real abfac(ncompdiff)
113      real,dimension(:,:),allocatable,save :: dij
114      real,save :: step
115
116
117! Initializations at first call
118      if (firstcall) then
119
120! Liste des especes qui sont diffuses si elles sont presentes
121!       ListeDiff=['co2','o','n2','ar','co','h2','h','d2','d','hd','o2','oh','ho2','h2o_vap','h2o2','o1d','n','no','no2']
122        ListeDiff(1)='co2'
123        ListeDiff(2)='o'
124        ListeDiff(3)='n2'
125        ListeDiff(4)='ar'
126        ListeDiff(5)='co'
127        ListeDiff(6)='h2'
128        ListeDiff(7)='h'
129        ListeDiff(8)='d2'
130        ListeDiff(9)='hd'
131        ListeDiff(10)='d'
132        ListeDiff(11)='o2'
133        ListeDiff(12)='h2o_vap'
134        ListeDiff(13)='o3'
135        ListeDiff(14)='n'
136        ListeDiff(15)='he'
137        ListeDiff(16)='hdo_vap'
138        i_h=1000
139        i_h2=1000
140        i_d=1000
141! On regarde quelles especes diffusables sont presentes
142
143        ncompdiff=0
144        indic_diff(1:nq)=0
145       
146        do nn=1,nq
147        do n=1,ListeDiffNb
148        if (ListeDiff(n) .eq. noms(nn)) then
149        indic_diff(nn)=1
150!       if (noms(nn) .eq. 'h') i_h=n
151!       if (noms(nn) .eq. 'h2') i_h2=n
152        endif
153
154        enddo
155        if (indic_diff(nn) .eq. 1) then
156        print*,'specie ', noms(nn), 'diffused in moldiff_red'
157        ncompdiff=ncompdiff+1
158        endif
159        enddo
160        print*,'number of diffused species:',ncompdiff
161        allocate(gcmind(ncompdiff))
162
163! Store gcm indexes in gcmind
164        n=0
165        do nn=1,nq
166        if (indic_diff(nn) .eq. 1) then
167        n=n+1
168        gcmind(n)=nn
169        if (noms(nn) .eq. 'h') i_h=n
170        if (noms(nn) .eq. 'h2') i_h2=n
171        if (noms(nn) .eq. 'd') i_d=n
172        endif
173        enddo
174
175!       print*,'gcmind',gcmind,i_h,i_h2
176
177! find vertical index above which diffusion is computed
178
179        do l=1,nlayer
180        if (pplay(1,l) .gt. Pdiff) then
181        il0=l
182        endif
183        enddo
184        il0=il0+1
185        print*,'vertical index for diffusion',il0,pplay(1,il0)
186
187        allocate(dij(ncompdiff,ncompdiff))
188
189        call moldiffcoeff_red(dij,indic_diff,gcmind,ncompdiff)
190        print*,'MOLDIFF  EXO'
191
192! allocatation des tableaux dependants du nombre d especes diffusees
193        allocate(qq(nlayer,ncompdiff))
194        allocate(qnew(nlayer,ncompdiff))
195        allocate(qint(nlayer,ncompdiff))
196        allocate(FacMass(nlayer,ncompdiff))
197        allocate(rhok(nlayer,ncompdiff))
198        allocate(rhokinit(nlayer,ncompdiff))
199
200        allocate(wi(ncompdiff))
201        allocate(wad(ncompdiff))
202        allocate(uthermal(ncompdiff))
203        allocate(lambdaexo(ncompdiff))
204        allocate(Hspecie(ncompdiff))
205        allocate(Mtot1(ncompdiff))
206        allocate(Mtot2(ncompdiff))
207        allocate(Mraf1(ncompdiff))
208        allocate(Mraf2(ncompdiff))
209
210        firstcall= .false.
211        step=1
212      endif ! of if (firstcall)
213
214!
215!ccccccccccccccccccccccccccccccccccccccccccccccccccccccc
216
217      masseU=1.660538782d-27
218      kbolt=1.3806504d-23
219      RgazP=8.314472   
220!      PI =3.141592653
221!      g=3.72d0
222      Rmars=3390000d0 ! Used to compute escape parameter no need to be more accurate
223        Grav=6.67d-11
224        Mmars=6.4d23
225        ij0=6000 ! For test
226        invsgmu=1d0/g/masseU
227       
228! Compute the wup(ig) for H and H2 using the balistic code from R Yelle
229
230        PhiEscH=0D0
231        PhiEscH2=0D0
232        PhiEscD=0D0
233
234      do ig=1,ngrid
235        pp=dble(pplay(ig,:))
236
237! Update the temperature
238
239!       CALL TMNEW(pt(ig,:),pdt(ig,:),pdtconduc(ig,:),pdteuv(ig,:)      &
240!     &  ,tt,ptimestep,nlayer,ig)
241        do l=1,nlayer
242          tt(l)=pt(ig,l)*1D0+(pdt(ig,l)*dble(ptimestep)+                &
243                              pdtconduc(ig,l)*dble(ptimestep)+          &
244                              pdteuv(ig,l)*dble(ptimestep))
245          ! to cach Nans...
246          if (tt(l).ne.tt(l)) then
247            print*,'Err TMNEW',ig,l,tt(l),pt(ig,l),                     &
248              pdt(ig,l),pdtconduc(ig,l),pdteuv(ig,l),dble(ptimestep)
249          endif
250        enddo ! of do l=1,nlayer
251
252! Update the mass mixing ratios modified by other processes
253
254!       CALL QMNEW(pq(ig,:,:),pdq(ig,:,:),qq,ptimestep,nlayer,        &
255!     &  ncompdiff,gcmind,ig)
256        do iq=1,ncompdiff
257         do l=1,nlayer
258          qq(l,iq)=pq(ig,l,gcmind(iq))*1D0+(                            &
259                           pdq(ig,l,gcmind(iq))*dble(ptimestep))
260          qq(l,iq)=max(qq(l,iq),1d-30)
261         enddo ! of do l=1,nlayer
262        enddo ! of do iq=1,ncompdiff
263
264! Compute the Pressure scale height
265
266        CALL HSCALE(pp,hp,nlayer)
267
268! Compute the atmospheric mass (in Dalton)
269
270        CALL MMOY(massemoy,mmol,qq,gcmind,nlayer,ncompdiff)
271
272! Compute the vertical gradient of atmospheric mass
273
274        CALL DMMOY(massemoy,hp,dmmeandz,nlayer)
275
276! Compute the altitude of each layer
277
278        CALL ZVERT(pp,tt,massemoy,zz,nlayer,ig)
279
280! Compute the total mass density (kg/m3)
281       
282        CALL RHOTOT(pp,tt,massemoy,qq,RHOT,RHOK,nlayer,ncompdiff)
283        RHOKINIT=RHOK
284
285! Compute total mass of each specie before new grid
286! For conservation tests
287! The conservation is not always fulfilled especially
288! for species very far from diffusion equilibrium "photochemical species"
289! e.g. OH, O(1D)
290
291        Mtot1(1:ncompdiff)=0d0
292
293        do l=il0,nlayer
294        do nn=1,ncompdiff
295        Mtot1(nn)=Mtot1(nn)+1d0/g*qq(l,nn)*                             &
296     &  (dble(pplev(ig,l))-dble(pplev(ig,l+1)))
297        enddo
298        enddo
299
300        Zmin=zz(il0)
301        Zmax=zz(nlayer)
302
303
304! If Zmax > 4000 km there is a problem / stop
305
306        if (Zmax .gt. 4000000.) then
307        Print*,'Zmax too high',ig,zmax,zmin
308        do l=1,nlayer
309        print*,'old',zz(l),pt(ig,l),pdteuv(ig,l),pdq(ig,l,:)
310        print*,'l',l,rhot(l),tt(l),pp(l),massemoy(l),qq(l,:)
311        enddo
312        stop
313        endif
314
315! The number of diffusion layers is variable
316! and fixed by the resolution (dzres) specified in diffusion.h
317! I fix a minimum  number of layers 40
318
319        nlraf=int((Zmax-Zmin)/1000./dzres)+1
320        if (nlraf .le. 40)  nlraf=40
321
322!       if (nlraf .ge. 200) print*,ig,nlraf,Zmin,Zmax
323       
324         ! allocate arrays:
325         allocate(Praf(nlraf),Traf(nlraf),Rraf(nlraf),Mraf(nlraf))
326         allocate(Nraf(nlraf),Draf(nlraf),Hraf(nlraf),Wraf(nlraf))
327         allocate(Zraf(nlraf),Tdiffraf(nlraf))
328         allocate(Prafold(nlraf),Mrafold(nlraf))
329         allocate(Qraf(nlraf,ncompdiff),Rrafk(nlraf,ncompdiff),Nrafk(nlraf,ncompdiff))
330         allocate(Rrafkold(nlraf,ncompdiff))
331         allocate(Drafmol(nlraf,ncompdiff),Hrafmol(nlraf,ncompdiff))
332         allocate(Wrafmol(nlraf,ncompdiff),Tdiffrafmol(nlraf,ncompdiff))
333         allocate(Atri(nlraf),Btri(nlraf),Ctri(nlraf),Dtri(nlraf),Xtri(nlraf))
334         allocate(Tad(nlraf),Dad(nlraf),Zad(nlraf),rhoad(nlraf))
335         allocate(alpha(nlraf),beta(nlraf),gama(nlraf),delta(nlraf),eps(nlraf))
336
337! before beginning, I use a better vertical resolution above il0,
338! altitude grid reinterpolation
339! The diffusion is solved on an altitude grid with constant step dz
340
341        CALL UPPER_RESOL(pp,tt,zz,massemoy,RHOT,RHOK,                   &
342     &  qq,mmol,gcmind,Praf,Traf,Qraf,Mraf,Zraf,                        &
343     &  Nraf,Nrafk,Rraf,Rrafk,il0,nlraf,ncompdiff,nlayer,ig)
344
345        Prafold=Praf
346        Rrafkold=Rrafk
347        Mrafold=Mraf
348       
349! We reddo interpolation of the gcm grid to estimate mass loss due to interpolation processes.
350
351        CALL GCMGRID_P(Zraf,Praf,Qraf,Traf,Nrafk,Rrafk,qq,qint,tt,tint  &
352     &    ,pp,mmol,gcmind,nlraf,ncompdiff,nlayer,ig)
353
354! We compute the mass correction factor of each specie at each pressure level
355
356        CALL CORRMASS(qq,qint,FacMass,nlayer,ncompdiff)
357
358! Altitude step
359
360        Dzraf=Zraf(2)-Zraf(1)
361
362! Total mass computed on the altitude grid
363
364        Mraf1(1:ncompdiff)=0d0
365        do nn=1,ncompdiff
366        do l=1,nlraf
367        Mraf1(nn)=Mraf1(nn)+Rrafk(l,nn)*Dzraf
368        enddo
369        enddo
370
371! Reupdate values for mass conservation
372
373!       do l=1,nlraf
374!       print*,'test',l,Nraf(l),Praf(l)
375!       do nn=1,ncompdiff
376!       Rrafk(l,nn)=RrafK(l,nn)*Mtot1(nn)/Mraf1(nn)
377!       enddo
378!       Rraf(l)=sum(Rrafk(l,:))
379!       do nn=1,ncompdiff
380!       Qraf(l,nn)=Rrafk(l,nn)/Rraf(l)
381!       Nrafk(l,nn)=Rrafk(l,nn)/dble(mmol(gcmind(nn)))/masseU
382!       enddo
383!       Mraf(l)=1d0/sum(Qraf(l,:)/dble(mmol(gcmind(:))))
384!       Nraf(l)=Rraf(l)/Mraf(l)/masseU
385!       Praf(l)=kbolt*Traf(l)*Nraf(l)
386!       print*,'test',l,Nraf(l),Praf(l)
387!       enddo
388
389!       do l=1,nlayer
390!       print*,'l',l,zz(l),pp(l),tt(l),sum(qq(l,:)),massemoy(l)
391!       enddo
392
393! The diffusion is computed above il0 computed from Pdiff in diffusion.h
394! No change below il0
395       
396        do l=1,nlayer
397         qnew(l,:)=qq(l,:) ! No effet below il0
398       enddo
399
400! all species treated independently
401
402! Upper boundary velocity
403! Jeans escape for H and H2
404! Comparison with and without escape still to be done
405! No escape for other species
406
407
408        do nn=1,ncompdiff
409        Uthermal(nn)=sqrt(2d0*kbolt*Traf(nlraf)/masseU/                 &
410     &  dble(mmol(gcmind(nn))))
411        Lambdaexo(nn)=masseU*dble(mmol(gcmind(nn)))*Grav*Mmars/         &
412     &  (Rmars+Zraf(nlraf))/kbolt/Traf(nlraf)
413        wi(nn)=0D0
414        if (nn .eq. i_h .or. nn .eq. i_h2 .or. nn .eq. i_d) then
415        wi(nn)=Uthermal(nn)/2d0/sqrt(PI)*exp(-Lambdaexo(nn))*           &
416     &  (Lambdaexo(nn)+1d0)
417        endif
418        enddo
419
420!       print*,'wi',wi(i_h),wi(i_h2),wi(i_d),Uthermal,Lambdaexo,mmol(gcmind(:))
421!       print*,'wi',wi
422
423! Compute time step for diffusion
424
425! Loop on species
426
427        T0=Traf(nlraf)
428        rho0=1d0
429
430        do nn=1,ncompdiff
431        masse=dble(mmol(gcmind(nn)))
432! DIffusion coefficient
433        CALL DCOEFF(nn,dij,Praf,Traf,Nraf,Nrafk,Draf,nlraf,ncompdiff)
434        Drafmol(:,nn)=Draf(:)
435! Scale height of the density of the specie
436        CALL HSCALEREAL(nn,Nrafk,Dzraf,Hraf,nlraf,ncompdiff)
437        Hrafmol(:,nn)=Hraf(:)
438!       Hspecie(nn)=kbolt*T0/masse*invsgmu
439! Computation of the diffusion vertical velocity of the specie
440        CALL VELVERT(nn,Traf,Hraf,Draf,Dzraf,masse,Wraf,nlraf)
441        Wrafmol(:,nn)=Wraf(:)
442! Computation of the diffusion time
443        CALL TIMEDIFF(nn,Hraf,Wraf,Tdiffraf,nlraf)
444        Tdiffrafmol(:,nn)=Tdiffraf
445        enddo
446! We use a lower time step
447        Tdiff=minval(Tdiffrafmol)
448        Tdiff=minval(Tdiffrafmol(nlraf,:))*Mraf(nlraf)
449! Some problems when H is dominant
450! The time step is chosen function of atmospheric mass at higher level
451! It is not very nice
452
453!       if (ig .eq. ij0) then
454!       print*,'test',ig,tdiff,tdiffmin,minloc(Tdiffrafmol),minloc(Tdiffrafmol(nlraf,:))
455!       endif
456        if (tdiff .lt. tdiffmin*Mraf(nlraf)) tdiff=tdiffmin*Mraf(nlraf)
457
458! Number of time step
459        ntime=anint(dble(ptimestep)/tdiff)
460!       print*,'ptime',ig,ptimestep,tdiff,ntime,tdiffmin,Mraf(nlraf)
461! Adimensionned temperature
462
463        do l=1,nlraf
464        Tad(l)=Traf(l)/T0
465        enddo
466
467        do istep=1,ntime
468        do nn=1,ncompdiff
469
470        Draf(1:nlraf)=Drafmol(1:nlraf,nn)
471
472! Parameters to adimension the problem
473
474        H0=kbolt*T0/dble(mmol(gcmind(nn)))*invsgmu
475        D0=Draf(nlraf)
476        Time0=H0*H0/D0
477        Time=Tdiff/Time0
478
479! Adimensions
480
481        do l=1,nlraf
482        Dad(l)=Draf(l)/D0
483        Zad(l)=Zraf(l)/H0
484        enddo
485        Wad(nn)=wi(nn)*Time0/H0
486        DZ=Zad(2)-Zad(1)
487        FacEsc=exp(-wad(nn)*DZ/Dad(nlraf-1))
488
489        do l=1,nlraf
490        RhoAd(l)=Rrafk(l,nn)/rho0
491        enddo
492
493! Compute intermediary coefficients
494
495        CALL DIFFPARAM(Tad,Dad,DZ,alpha,beta,gama,delta,eps,nlraf)
496
497! First way to include escape need to be validated
498! Compute escape factor (exp(-ueff*dz/D(nz)))
499
500! Compute matrix coefficients
501
502        CALL MATCOEFF(alpha,beta,gama,delta,eps,Dad,rhoAd,FacEsc,       &
503     &   dz,time,Atri,Btri,Ctri,Dtri,nlraf)
504
505        Xtri(:)=0D0
506
507! SOLVE SYSTEM
508
509        CALL tridag(Atri,Btri,Ctri,Dtri,Xtri,nlraf)
510
511!       Xtri=rhoAd
512
513!       if (ig .eq. ij0 .and. (nn .eq. 1 .or. nn .eq. 5 .or. nn .eq. 6 .or. nn .eq. 16)) then
514!       do l=1,nlraf
515!       if (Xtri(l) .lt. 0.) then
516!       print*,'l',l,rhoAd(l)*rho0,Xtri(l)*rho0,nn,Tad(l),Zad(l),Dad(l)
517!       stop
518!       endif
519!       enddo
520!       endif
521
522! Check mass conservation of speci
523
524!       CALL CheckMass(rhoAd,Xtri,nlraf,nn)
525
526! Update mass density of the species
527
528        do l=1,nlraf
529        Rrafk(l,nn)=rho0*Xtri(l)
530        if (Rrafk(l,nn) .ne. Rrafk(l,nn) .or.                           &
531     &  Rrafk(l,nn) .lt. 0 .and. nn .eq. 16) then
532
533! Test if n(CO2) < 0 skip diffusion (should never happen)
534
535        print*,'pb moldiff',istep,ig,l,nn,Rrafk(l,nn),tdiff,&
536     &  rho0*Rhoad(l),Zmin,Zmax,ntime
537        print*,'Atri',Atri
538        print*,'Btri',Btri
539        print*,'Ctri',Ctri
540        print*,'Dtri',Dtri
541        print*,'Xtri',Xtri
542        print*,'alpha',alpha
543        print*,'beta',beta
544        print*,'gama',gama
545        print*,'delta',delta
546        print*,'eps',eps
547        print*,'Dad',Dad
548        print*,'Temp',Traf
549        print*,'alt',Zraf
550        print*,'Mraf',Mraf
551        stop
552!       pdqdiff(1:ngrid,1:nlayer,1:nq)=0.
553!       return
554!       Rrafk(l,nn)=1D-30*Rraf(l)
555        Rrafk(l,nn)=rho0*Rhoad(l)
556
557        endif
558
559        enddo
560
561        enddo ! END Species Loop
562
563! Update total mass density
564
565        do l=1,nlraf
566        Rraf(l)=sum(Rrafk(l,:))
567        enddo
568
569! Compute new mass average at each altitude level and new pressure
570
571        do l=1,nlraf
572        do nn=1,ncompdiff
573        Qraf(l,nn)=Rrafk(l,nn)/Rraf(l)
574        Nrafk(l,nn)=Rrafk(l,nn)/dble(mmol(gcmind(nn)))/masseU
575        enddo
576        Mraf(l)=1d0/sum(Qraf(l,:)/dble(mmol(gcmind(:))))
577        Nraf(l)=Rraf(l)/Mraf(l)/masseU
578        Praf(l)=Nraf(l)*kbolt*Traf(l)
579        enddo
580
581        enddo ! END time Loop
582
583! Compute the total mass of each species to check mass conservation     
584
585        Mraf2(1:ncompdiff)=0d0
586        do nn=1,ncompdiff
587        do l=1,nlraf
588        Mraf2(nn)=Mraf2(nn)+Rrafk(l,nn)*Dzraf
589        enddo
590        enddo
591
592!        print*,'Mraf',Mraf2
593
594! Reinterpolate values on the GCM pressure levels
595
596        CALL GCMGRID_P2(Zraf,Praf,Qraf,Traf,Nrafk,Rrafk,qq,qnew,tt,tnew,&
597     &   pp,mmol,gcmind,nlraf,ncompdiff,nlayer,FacMass,ig)
598
599        CALL RHOTOT(pp,tt,massemoy,qnew,RHOT,RHOK,nlayer,ncompdiff)
600
601! Update total escape flux of H and H2 (if q was really qnew, but not forget we will output
602! the trend only at the end
603
604!       if (i_h .ne. 1000) PhiEscH=PhiEscH+wi(i_h)*Nrafk(nlraf,i_h)*cell_area(ig) ! in s-1
605        if (i_h .ne. 1000) PhiauxH(ig)=wi(i_h)*Nrafk(nlraf,i_h)*cell_area(ig) ! in s-1
606!       if (i_h2 .ne. 1000) PhiEscH2=PhiEscH2+wi(i_h2)*Nrafk(nlraf,i_h2)*cell_area(ig) ! in s-1 (U in m/s, aire in m2, Nrafk in m-3)
607        if (i_h2 .ne. 1000) PhiauxH2(ig)=wi(i_h2)*Nrafk(nlraf,i_h2)*cell_area(ig)
608!        if (i_d .ne. 1000) PhiEscD=PhiEscD+wi(i_d)*Nrafk(nlraf,i_d)*cell_area(ig)
609        if (i_d .ne. 1000) PhiauxD(ig)=wi(i_d)*Nrafk(nlraf,i_d)*cell_area(ig)
610!       print*,'test',ig,wi(i_h),wi(i_h2),Nrafk(nlraf,i_h),Nrafk(nlraf,i_h2),Nrafk(nlraf,i_d),cell_area(ig),PhiEscH,PhiEscH2,i_h,i_h2,i_d,PhiEscD
611!       stop
612
613
614        if (ig .eq. ij0) then
615        do l=il0,nlayer
616        write(*,'(i2,1x,19(e12.4,1x))') l,zz(l),tt(l),RHOK(l,1)/sum(RHOK(l,:)),RHOKINIT(l,1)/sum(RHOKINIT(l,:)),&
617     &  RHOK(l,2)/sum(RHOK(l,:)),RHOKINIT(l,2)/sum(RHOKINIT(l,:)),&
618     &  RHOK(l,6)/sum(RHOK(l,:)),RHOKINIT(l,6)/sum(RHOKINIT(l,:)),&
619     &  RHOK(l,5)/sum(RHOK(l,:)),RHOKINIT(l,5)/sum(RHOKINIT(l,:)),&
620     &  RHOK(l,7)/sum(RHOK(l,:)),RHOKINIT(l,7)/sum(RHOKINIT(l,:))
621        enddo
622        endif
623
624! Compute total mass of each specie on the GCM vertical grid
625
626        Mtot2(1:ncompdiff)=0d0
627
628        do l=il0,nlayer
629        do nn=1,ncompdiff
630        Mtot2(nn)=Mtot2(nn)+1d0/g*qnew(l,nn)*                           &
631     &  (dble(pplev(ig,l))-dble(pplev(ig,l+1)))
632        enddo
633        enddo
634
635! Check mass  conservation of each specie on column
636
637!       do nn=1,ncompdiff
638!       CALL CheckMass2(qq,qnew,pplev(ig,:),il0,nlayer,nn,ncompdiff)
639!       enddo
640
641! Compute the diffusion trends du to diffusion
642
643        do l=1,nlayer
644        do nn=1,ncompdiff
645        pdqdiff(ig,l,gcmind(nn))=(qnew(l,nn)-qq(l,nn))/ptimestep
646        enddo
647        enddo
648
649! deallocation des tableaux
650
651        deallocate(Praf,Traf,Rraf,Mraf)
652        deallocate(Nraf,Draf,Hraf,Wraf)
653        deallocate(Zraf,Tdiffraf)
654        deallocate(Prafold,Mrafold)
655        deallocate(Qraf,Rrafk,Nrafk)
656        deallocate(Rrafkold)
657        deallocate(Drafmol,Hrafmol)
658        deallocate(Wrafmol,Tdiffrafmol)
659        deallocate(Atri,Btri,Ctri,Dtri,Xtri)
660        deallocate(Tad,Dad,Zad,rhoad)
661        deallocate(alpha,beta,gama,delta,eps)
662
663
664       enddo  ! ig loop         
665!       print*,'Escape flux H, H2,D (s-1)',PhiEscH,PhiEscH2,PhiEscD
666       if (i_h.ne.1000) call planetwide_sumval(PhiauxH,PhiEscH)
667       if (i_h2.ne.1000) call planetwide_sumval(PhiauxH2,PhiEscH2)
668       if (i_d.ne.1000) call planetwide_sumval(PhiauxD,PhiEscD)
669
670      return
671      end
672
673!    ********************************************************************
674!    ********************************************************************
675!    ********************************************************************
676 
677!     JYC subtroutine solving MX = Y where M is defined as a block tridiagonal
678!       matrix (Thomas algorithm), tested on a example
679
680      subroutine tridagbloc(M,F,X,n1,n2)
681      parameter (nmax=100)
682      real*8 M(n1*n2,n1*n2),F(n1*n2),X(n1*n2)
683      real*8 A(n1,n1,n2),B(n1,n1,n2),C(n1,n1,n2),D(n1,n2)
684      real*8 at(n1,n1),bt(n1,n1),ct(n1,n1),dt(n1),gamt(n1,n1),y(n1,n1)
685      real*8 alf(n1,n1),gam(n1,n1,n2),alfinv(n1,n1)
686      real*8 uvec(n1,n2),uvect(n1),vvect(n1),xt(n1)
687      real*8 indx(n1)
688      real*8 rhu
689      integer n1,n2
690      integer i,p,q   
691
692        X(:)=0.
693! Define the bloc matrix A,B,C and the vector D
694      A(1:n1,1:n1,1)=M(1:n1,1:n1)
695      C(1:n1,1:n1,1)=M(1:n1,n1+1:2*n1)
696      D(1:n1,1)=F(1:n1)
697
698      do i=2,n2-1
699      A(1:n1,1:n1,i)=M((i-1)*n1+1:i*n1,(i-1)*n1+1:i*n1)
700      B(1:n1,1:n1,i)=M((i-1)*n1+1:i*n1,(i-2)*n1+1:(i-1)*n1)
701      C(1:n1,1:n1,i)=M((i-1)*n1+1:i*n1,i*n1+1:(i+1)*n1)
702      D(1:n1,i)=F((i-1)*n1+1:i*n1)
703      enddo
704      A(1:n1,1:n1,n2)=M((n2-1)*n1+1:n2*n1,(n2-1)*n1+1:n2*n1)
705      B(1:n1,1:n1,n2)=M((n2-1)*n1+1:n2*n1,(n2-2)*n1+1:(n2-1)*n1)
706      D(1:n1,n2)=F((n2-1)*n1+1:n2*n1)
707
708! Initialization
709        y(:,:)=0.
710        do i=1,n1
711        y(i,i)=1.
712        enddo
713       
714        at(:,:)=A(:,:,1)
715        ct(:,:)=C(:,:,1)
716        dt(:)=D(:,1)
717        call ludcmp(at,n1,n1,indx,rhu,ierr)
718        do p=1,n1
719         call lubksb(at,n1,n1,indx,y(1,p))
720         do q=1,n1
721         alfinv(q,p)=y(q,p)
722         enddo
723        enddo
724      gamt=matmul(alfinv,ct)
725      gam(:,:,1)=gamt(:,:)     
726      uvect=matmul(alfinv,dt)
727      uvec(:,1)=uvect
728
729      do i=2,n2-1
730        y(:,:)=0.
731       do j=1,n1
732         y(j,j)=1.
733       enddo 
734      bt(:,:)=B(:,:,i)
735      at(:,:)=A(:,:,i)-matmul(bt,gamt)
736      ct(:,:)=C(:,:,i)
737      dt(:)=D(:,i)
738      call ludcmp(at,n1,n1,indx,rhu,ierr)
739        do p=1,n1
740         call lubksb(at,n1,n1,indx,y(1,p))
741         do q=1,n1
742         alfinv(q,p)=y(q,p)
743         enddo
744        enddo
745      gamt=matmul(alfinv,ct)
746      gam(:,:,i)=gamt
747      vvect=dt-matmul(bt,uvect)
748      uvect=matmul(alfinv,vvect)
749      uvec(:,i)=uvect
750      enddo
751      bt=B(:,:,n2)
752      dt=D(:,n2)
753      at=A(:,:,n2)-matmul(bt,gamt)
754      vvect=dt-matmul(bt,uvect)
755       y(:,:)=0.
756       do j=1,n1
757         y(j,j)=1.
758       enddo
759       call ludcmp(at,n1,n1,indx,rhu,ierr)
760        do p=1,n1
761         call lubksb(at,n1,n1,indx,y(1,p))
762         do q=1,n1
763         alfinv(q,p)=y(q,p)
764         enddo
765        enddo
766      xt=matmul(alfinv,vvect)
767      X((n2-1)*n1+1 :n1*n2)=xt
768      do i=n2-1,1,-1
769      gamt=gam(:,:,i)
770      xt=X(i*n1+1:n1*n2)
771      uvect=uvec(:,i)
772      vvect=matmul(gamt,xt)
773      X((i-1)*n1+1:i*n1)=uvect-vvect
774      enddo
775
776      end
777
778      subroutine tridag(a,b,c,r,u,n)
779!      parameter (nmax=4000)
780!      dimension gam(nmax),a(n),b(n),c(n),r(n),u(n)
781      real*8 gam(n),a(n),b(n),c(n),r(n),u(n)
782      if(b(1).eq.0.)then
783        stop 'tridag: error: b(1)=0 !!! '
784      endif
785      bet=b(1)
786      u(1)=r(1)/bet
787      do 11 j=2,n
788        gam(j)=c(j-1)/bet
789        bet=b(j)-a(j)*gam(j)
790        if(bet.eq.0.) then
791          stop 'tridag: error: bet=0 !!! '
792        endif
793        u(j)=(r(j)-a(j)*u(j-1))/bet
79411    continue
795      do 12 j=n-1,1,-1
796        u(j)=u(j)-gam(j+1)*u(j+1)
79712    continue
798      return
799      end
800
801!    ********************************************************************
802!    ********************************************************************
803!    ********************************************************************
804
805      SUBROUTINE LUBKSB(A,N,NP,INDX,B)
806
807      implicit none
808
809      integer i,j,n,np,ii,ll
810      real*8 sum
811      real*8 a(np,np),indx(np),b(np)
812
813!      DIMENSION A(NP,NP),INDX(N),B(N)
814      II=0
815      DO 12 I=1,N
816        LL=INDX(I)
817        SUM=B(LL)
818        B(LL)=B(I)
819        IF (II.NE.0)THEN
820          DO 11 J=II,I-1
821            SUM=SUM-A(I,J)*B(J)
82211        CONTINUE
823        ELSE IF (SUM.NE.0.) THEN
824          II=I
825        ENDIF
826        B(I)=SUM
82712    CONTINUE
828      DO 14 I=N,1,-1
829        SUM=B(I)
830        IF(I.LT.N)THEN
831          DO 13 J=I+1,N
832            SUM=SUM-A(I,J)*B(J)
83313        CONTINUE
834        ENDIF
835        B(I)=SUM/A(I,I)
83614    CONTINUE
837      RETURN
838      END
839
840!    ********************************************************************
841!    ********************************************************************
842!    ********************************************************************
843
844      SUBROUTINE LUDCMP(A,N,NP,INDX,D,ierr)
845
846      implicit none
847
848      integer n,np,nmax,i,j,k,imax
849      real*8 d,tiny,aamax
850      real*8 a(np,np),indx(np)
851      integer ierr  ! error =0 if OK, =1 if problem
852
853      PARAMETER (NMAX=100,TINY=1.0E-20)
854!      DIMENSION A(NP,NP),INDX(N),VV(NMAX)
855      real*8 sum,vv(nmax),dum
856
857      D=1.
858      DO 12 I=1,N
859        AAMAX=0.
860        DO 11 J=1,N
861          IF (ABS(A(I,J)).GT.AAMAX) AAMAX=ABS(A(I,J))
86211      CONTINUE
863        IF (AAMAX.EQ.0.) then
864            write(*,*) 'In moldiff: Problem in LUDCMP with matrix A'
865            write(*,*) 'Singular matrix ?'
866            write(*,*) 'Matrix A = ', A
867!            stop
868!           TO DEBUG :
869            ierr =1
870            return
871!           stop
872        END IF
873
874        VV(I)=1./AAMAX
87512    CONTINUE
876      DO 19 J=1,N
877        IF (J.GT.1) THEN
878          DO 14 I=1,J-1
879            SUM=A(I,J)
880            IF (I.GT.1)THEN
881              DO 13 K=1,I-1
882                SUM=SUM-A(I,K)*A(K,J)
88313            CONTINUE
884              A(I,J)=SUM
885            ENDIF
88614        CONTINUE
887        ENDIF
888        AAMAX=0.
889        DO 16 I=J,N
890          SUM=A(I,J)
891          IF (J.GT.1)THEN
892            DO 15 K=1,J-1
893              SUM=SUM-A(I,K)*A(K,J)
89415          CONTINUE
895            A(I,J)=SUM
896          ENDIF
897          DUM=VV(I)*ABS(SUM)
898          IF (DUM.GE.AAMAX) THEN
899            IMAX=I
900            AAMAX=DUM
901          ENDIF
90216      CONTINUE
903        IF (J.NE.IMAX)THEN
904          DO 17 K=1,N
905            DUM=A(IMAX,K)
906            A(IMAX,K)=A(J,K)
907            A(J,K)=DUM
90817        CONTINUE
909          D=-D
910          VV(IMAX)=VV(J)
911        ENDIF
912        INDX(J)=IMAX
913        IF(J.NE.N)THEN
914          IF(A(J,J).EQ.0.)A(J,J)=TINY
915          DUM=1./A(J,J)
916          DO 18 I=J+1,N
917            A(I,J)=A(I,J)*DUM
91818        CONTINUE
919        ENDIF
92019    CONTINUE
921      IF(A(N,N).EQ.0.)A(N,N)=TINY
922      ierr =0
923      RETURN
924      END
925
926        SUBROUTINE TMNEW(T1,DT1,DT2,DT3,T2,dtime,nl,ig)
927        IMPLICIT NONE
928
929        INTEGER,INTENT(IN) :: nl,ig
930        REAL,INTENT(IN),DIMENSION(nl) :: T1,DT1,DT2,DT3
931        REAL*8,INTENT(OUT),DIMENSION(nl) :: T2
932        REAL,INTENT(IN) :: dtime
933        INTEGER :: l
934        DO l=1,nl
935        T2(l)=T1(l)*1D0+(DT1(l)*dble(dtime)+                            &
936     &  DT2(l)*dble(dtime)+                                             &
937     &  DT3(l)*dble(dtime))*1D0
938        if (T2(l) .ne. T2(l)) then
939        print*,'Err TMNEW',ig,l,T2(l),T1(l),dT1(l),DT2(l),              &
940     &  DT3(l),dtime,dble(dtime)
941        endif
942
943        ENDDO
944        END
945
946        SUBROUTINE QMNEW(Q1,DQ,Q2,dtime,nl,nq,gc,ig)
947        use tracer_mod, only: nqmx
948        IMPLICIT NONE
949
950        INTEGER,INTENT(IN) :: nl,nq
951        INTEGER,INTENT(IN) :: ig
952        INTEGER,INTENT(IN),dimension(nq) :: gc
953        REAL,INTENT(IN),DIMENSION(nl,nqmx) :: Q1,DQ
954        REAL*8,INTENT(OUT),DIMENSION(nl,nq) :: Q2
955        REAL,INTENT(IN) :: dtime
956        INTEGER :: l,iq
957        DO l=1,nl
958        DO iq=1,nq
959        Q2(l,iq)=Q1(l,gc(iq))*1D0+(DQ(l,gc(iq))*dble(dtime))*1D0
960        Q2(l,iq)=max(Q2(l,iq),1d-30)
961        ENDDO
962        ENDDO
963        END
964
965        SUBROUTINE HSCALE(p,hp,nl)
966        IMPLICIT NONE
967
968        INTEGER :: nl
969        INTEGER :: l
970        REAL*8,dimension(nl) :: P
971        REAL*8,DIMENSION(nl) :: Hp
972       
973        hp(1)=-log(P(2)/P(1))
974        hp(nl)=-log(P(nl)/P(nl-1))
975
976        DO l=2,nl-1
977        hp(l)=-log(P(l+1)/P(l-1))
978        ENDDO
979        END
980
981        SUBROUTINE MMOY(massemoy,mmol,qq,gc,nl,nq)
982        use tracer_mod, only: nqmx
983        IMPLICIT NONE
984
985        INTEGER :: nl,nq,l
986        INTEGER,dimension(nq) :: gc
987        REAL*8,DIMENSION(nl,nq) :: qq
988        REAL*8,DIMENSION(nl) :: massemoy
989        REAL,DIMENSION(nqmx) :: MMOL
990
991
992        do l=1,nl
993        massemoy(l)=1D0/sum(qq(l,:)/dble(mmol(gc(:))))
994        enddo
995
996        END
997
998        SUBROUTINE DMMOY(M,H,DM,nl)
999        IMPLICIT NONE   
1000        INTEGER :: nl,l
1001        REAL*8,DIMENSION(nl) :: M,H,DM
1002
1003        DM(1)=(-3D0*M(1)+4D0*M(2)-M(3))/2D0/H(1)
1004        DM(nl)=(3D0*M(nl)-4D0*M(nl-1)+M(nl-2))/2D0/H(nl)
1005
1006        do l=2,nl-1
1007        DM(l)=(M(l+1)-M(l-1))/H(l)
1008        enddo
1009
1010        END
1011
1012        SUBROUTINE ZVERT(P,T,M,Z,nl,ig)
1013        IMPLICIT NONE
1014        INTEGER :: nl,l,ig
1015        REAL*8,dimension(nl) :: P,T,M,Z,H
1016        REAL*8 :: z0
1017        REAL*8 :: kbolt,masseU,Konst,g,Hpm
1018        masseU=1.660538782d-27
1019        kbolt=1.3806504d-23
1020        Konst=kbolt/masseU
1021        g=3.72D0
1022
1023        z0=0d0
1024        Z(1)=z0
1025        H(1)=Konst*T(1)/M(1)/g
1026
1027        do l=2,nl
1028        H(l)=Konst*T(l)/M(l)/g
1029        Hpm=H(l-1)
1030        Z(l)=z(l-1)-Hpm*log(P(l)/P(l-1))
1031        if (Z(l) .ne. Z(l)) then
1032        print*,'pb',l,ig
1033        print*,'P',P
1034        print*,'T',T
1035        print*,'M',M
1036        print*,'Z',Z
1037        print*,'Hpm',Hpm
1038        endif
1039        enddo
1040
1041        END
1042
1043        SUBROUTINE RHOTOT(P,T,M,qq,rhoN,rhoK,nl,nq)
1044        IMPLICIT NONE
1045
1046        REAL*8 :: kbolt,masseU,Konst
1047        INTEGER :: nl,nq,l,iq
1048        REAL*8,DIMENSION(nl) :: P,T,M,RHON
1049        REAL*8,DIMENSION(nl,nq) :: RHOK,qq
1050        masseU=1.660538782d-27
1051        kbolt=1.3806504d-23
1052        Konst=Kbolt/masseU
1053
1054        do l=1,nl
1055        RHON(l)=P(l)*M(l)/T(l)/Konst
1056        do iq=1,nq
1057        RHOK(l,iq)=qq(l,iq)*RHON(l)
1058        enddo
1059        enddo
1060
1061        END
1062
1063        SUBROUTINE UPPER_RESOL(P,T,Z,M,R,Rk,                            &
1064     & qq,mmol,gc,Praf,Traf,Qraf,Mraf,Zraf,                             &
1065     & Nraf,Nrafk,Rraf,Rrafk,il,nl,nq,nlx,ig)
1066        use tracer_mod, only: nqmx
1067        IMPLICIT NONE
1068       
1069        INTEGER :: nl,nq,il,l,i,iq,nlx,iz,ig
1070        INTEGER :: gc(nq)
1071        INTEGER,DIMENSION(1) :: indz
1072        REAL*8, DIMENSION(nlx) :: P,T,Z,M,R
1073        REAL*8, DIMENSION(nlx,nq) :: qq,Rk
1074        REAL*8, DIMENSION(nl) :: Praf,Traf,Mraf,Zraf,Nraf,Rraf
1075        REAL*8 :: kbolt,masseU,Konst,g
1076        REAL*8, DIMENSION(nl,nq) :: Qraf,Rrafk,Nrafk
1077        REAL*8 :: facZ,dZ,H
1078        REAL,DIMENSION(nqmx) :: mmol
1079        masseU=1.660538782d-27
1080        kbolt=1.3806504d-23
1081        Konst=Kbolt/masseU
1082        g=3.72d0
1083
1084
1085        Zraf(1)=z(il)
1086        Praf(1)=P(il)
1087        Traf(1)=T(il)
1088        Nraf(1)=Praf(1)/kbolt/Traf(1)
1089        do iq=1,nq
1090        Qraf(1,iq)=qq(il,iq)
1091        enddo
1092        Mraf(1)=1d0/sum(Qraf(1,:)/dble(mmol(gc(:))))
1093        Rraf(1)=Mraf(1)*masseU*Nraf(1)
1094        do iq=1,nq
1095        Rrafk(1,iq)=Rraf(1)*Qraf(1,iq)
1096        Nrafk(1,iq)=Rrafk(1,iq)/masseU/dble(mmol(gc(iq)))
1097        enddo
1098        Zraf(nl)=z(nlx)
1099
1100        do l=2,nl-1
1101        Zraf(l)=Zraf(1)+(Zraf(nl)-Zraf(1))/dble(nl-1)*dble((l-1))
1102        indz=maxloc(z,mask=z <= Zraf(l))
1103        iz=indz(1)
1104        if (iz .lt. 1 .or. iz .gt. nlx) then
1105        print*,'danger',iz,nl,Zraf(l),l,Zraf(1),Zraf(nl),z,P,T,ig
1106        stop
1107        endif
1108        dZ=Zraf(l)-Zraf(l-1)
1109!       dZ=Zraf(l)-z(iz)
1110        facz=(Zraf(l)-z(iz))/(z(iz+1)-z(iz))
1111        Traf(l)=T(iz)+(T(iz+1)-T(iz))*facz
1112        do iq=1,nq
1113!        Qraf(l,iq)=qq(iz,iq)+(qq(iz+1,iq)-qq(iz,iq))*facz
1114        Rrafk(l,iq)=Rk(iz,iq)+(Rk(iz+1,iq)-Rk(iz,iq))*facZ
1115        Rrafk(l,iq)=Rk(iz,iq)*(Rk(iz+1,iq)/Rk(iz,iq))**facZ
1116        enddo
1117!       Mraf(l)=1D0/(sum(qraf(l,:)/dble(mmol(gc(:)))))
1118        Rraf(l)=sum(Rrafk(l,:))
1119        do iq=1,nq
1120        Qraf(l,iq)=Rrafk(l,iq)/Rraf(l)
1121        enddo
1122        Mraf(l)=1D0/(sum(qraf(l,:)/dble(mmol(gc(:)))))
1123!       H=Konst*Traf(l)/Mraf(l)/g
1124!       H=Konst*T(iz)/M(iz)/g
1125!       Praf(l)=P(iz)*exp(-dZ/H)
1126!       Praf(l)=Praf(l-1)*exp(-dZ/H)
1127!       print*,'iz',l,iz,Praf(il-1)*exp(-dZ/H),z(iz),z(iz+1),H
1128        Nraf(l)=Rraf(l)/Mraf(l)/masseU
1129        Praf(l)=Nraf(l)*kbolt*Traf(l)
1130!       Rraf(l)=Nraf(l)*Mraf(l)*masseU
1131        do iq=1,nq
1132!       Rrafk(l,iq)=Rraf(l)*Qraf(l,iq)
1133        Nrafk(l,iq)=Rrafk(l,iq)/dble(mmol(gc(iq)))/masseU
1134        if (Nrafk(l,iq) .lt. 0. .or.                                    &
1135     &  Nrafk(l,iq) .ne. Nrafk(l,iq)) then
1136        print*,'pb interpolation',l,iq,Nrafk(l,iq),Rrafk(l,iq),         &
1137     &  Qraf(l,iq),Rk(iz,iq),Rk(iz+1,iq),facZ,Zraf(l),z(iz)
1138        stop
1139        endif
1140        enddo
1141        enddo
1142        Zraf(nl)=z(nlx)
1143        Traf(nl)=T(nlx)
1144        do iq=1,nq
1145        Rrafk(nl,iq)=Rk(nlx,iq)
1146        Qraf(nl,iq)=Rk(nlx,iq)/R(nlx)
1147        Nrafk(nl,iq)=Rk(nlx,iq)/dble(mmol(gc(iq)))/masseU
1148        enddo
1149        Nraf(nl)=sum(Nrafk(nl,:))
1150        Praf(nl)=Nraf(nl)*kbolt*Traf(nl)
1151        Mraf(nl)=1D0/sum(Qraf(nl,:)/dble(mmol(gc(:))))
1152        END
1153
1154        SUBROUTINE CORRMASS(qq,qint,FacMass,nl,nq)
1155        IMPLICIT NONE
1156        INTEGER :: nl,nq,l,nn
1157        REAL*8,DIMENSION(nl,nq) :: qq,qint,FacMass
1158
1159        do nn=1,nq
1160        do l=1,nl
1161        FacMass(l,nn)=qq(l,nn)/qint(l,nn)
1162        enddo
1163        enddo
1164
1165        END     
1166
1167
1168        SUBROUTINE DCOEFF(nn,Dij,P,T,N,Nk,D,nl,nq)
1169        IMPLICIT NONE
1170        REAL*8,DIMENSION(nl) :: N,T,D,P
1171        REAL*8,DIMENSION(nl,nq) :: Nk
1172        INTEGER :: nn,nl,nq,l,iq
1173        REAL,DIMENSION(nq,nq)  :: Dij
1174        REAL*8 :: interm,P0,T0,ptfac,dfac
1175
1176        P0=1D5
1177        T0=273d0
1178       
1179
1180        do l=1,nl
1181        ptfac=(P0/P(l))*(T(l)/T0)**1.75d0
1182        D(l)=0d0
1183        interm=0d0
1184        do iq=1,nq
1185        if (iq .ne. nn) then
1186        dfac=dble(dij(nn,iq))*ptfac
1187        interm=interm+Nk(l,iq)/N(l)/dfac
1188        endif
1189        enddo
1190        !Temporary: eliminate modification to include Wilke's formulation
1191        !back to the old scheme to check effect
1192        !D(l)=1d0/interm
1193        D(l)=(1D0-Nk(l,nn)/N(l))/interm
1194        enddo
1195        END
1196 
1197        SUBROUTINE HSCALEREAL(nn,Nk,Dz,H,nl,nq)
1198        IMPLICIT NONE
1199        INTEGER :: nn,nl,nq,l
1200        REAL*8,DIMENSION(nl) :: H
1201        REAL*8,DIMENSION(nl,nq) :: Nk
1202        REAL*8 :: Dz
1203       
1204        H(1)=(-3D0*Nk(1,nn)+4d0*NK(2,nn)-Nk(3,nn))/(2D0*DZ)/            &
1205     &  NK(1,nn)
1206
1207        H(1)=-1D0/H(1)
1208
1209        DO l=2,nl-1
1210        H(l)=(Nk(l+1,nn)-NK(l-1,nn))/(2D0*DZ)/NK(l,nn)
1211        H(l)=-1D0/H(l)
1212        ENDDO
1213
1214        H(nl)=(3D0*Nk(nl,nn)-4D0*Nk(nl-1,nn)+Nk(nl-2,nn))/(2D0*DZ)/     &
1215     &  Nk(nl,nn)
1216        H(nl)=-1D0/H(nl)
1217
1218!       do l=1,nl
1219!       if (abs(H(l)) .lt. 100.) then
1220!       print*,'H',l,H(l),Nk(l,nn),nn
1221!       endif
1222!       enddo
1223
1224        END
1225       
1226        SUBROUTINE VELVERT(nn,T,H,D,Dz,masse,W,nl)
1227        IMPLICIT NONE
1228        INTEGER :: l,nl,nn
1229        REAL*8,DIMENSION(nl) :: T,H,D,W,DT
1230        REAL*8 :: Dz,Hmol,masse
1231        REAL*8 :: kbolt,masseU,Konst,g
1232        masseU=1.660538782d-27
1233        kbolt=1.3806504d-23
1234        Konst=Kbolt/masseU
1235        g=3.72d0
1236
1237        DT(1)=1D0/T(1)*(-3D0*T(1)+4D0*T(2)-T(3))/(2D0*DZ)
1238        Hmol=Konst*T(1)/masse/g
1239        W(1)=-D(1)*(1D0/H(1)-1D0/Hmol-DT(1))
1240
1241        DO l=2,nl-1
1242        DT(l)=1D0/T(l)*(T(l+1)-T(l-1))/(2D0*DZ)
1243        Hmol=Konst*T(l)/masse/g
1244        W(l)=-D(l)*(1D0/H(l)-1D0/Hmol-DT(l))
1245        ENDDO
1246
1247        DT(nl)=1D0/T(nl)*(3d0*T(nl)-4D0*T(nl-1)+T(nl-2))/(2D0*DZ)
1248        Hmol=Konst*T(nl)/masse/g
1249        W(nl)=-D(nl)*(1D0/H(nl)-1D0/Hmol-DT(nl))
1250
1251!       do l=1,nl
1252!       print*,'W',W(l),D(l),H(l),DT(l)
1253!       enddo
1254
1255        END
1256
1257        SUBROUTINE TIMEDIFF(nn,H,W,TIME,nl)
1258        IMPLICIT NONE
1259        INTEGER :: nn,nl,l
1260        REAL*8,DIMENSION(nl) :: W,H,TIME
1261
1262        DO l=1,nl
1263        TIME(l)=abs(H(l)/W(l))
1264        if (TIME(l) .lt. 1.D-4) then
1265!       print*,'low tdiff',H(l),W(l),nn,l
1266        endif
1267        ENDDO
1268
1269        END
1270
1271
1272        SUBROUTINE DIFFPARAM(T,D,dz,alpha,beta,gama,delta,eps,nl)
1273        IMPLICIT NONE
1274        INTEGER :: nl,l
1275        REAL*8,DIMENSION(nl) :: T,D
1276        REAL*8 :: DZ,DZinv
1277        REAL*8,DIMENSION(nl) :: alpha,beta,gama,delta,eps
1278
1279! Compute alpha,beta and delta
1280! lower altitude values
1281        DZinv=1d0/(2D0*DZ)
1282
1283        beta(1)=1d0/T(1)
1284        alpha(1)=beta(1)*(-3D0*T(1)+4D0*T(2)-T(3))*Dzinv
1285        delta(1)=(-3D0*D(1)+4D0*D(2)-D(3))*Dzinv
1286
1287        beta(2)=1d0/T(2)
1288        alpha(2)=beta(2)*(T(3)-T(1))*Dzinv
1289        delta(2)=(D(3)-D(1))*Dzinv
1290
1291!       do l=2,nl-1
1292!       beta(l)=1D0/T(l)
1293!       alpha(l)=beta(l)*(T(l+1)-T(l-1))*Dzinv
1294!       delta(l)=(D(l+1)-D(l-1))*Dzinv
1295!       end do   
1296
1297        do l=3,nl-1
1298        beta(l)=1D0/T(l)
1299        alpha(l)=beta(l)*(T(l+1)-T(l-1))*Dzinv
1300        delta(l)=(D(l+1)-D(l-1))*Dzinv
1301        gama(l-1)=(beta(l)-beta(l-2))*Dzinv
1302        eps(l-1)=(alpha(l)-alpha(l-2))*Dzinv
1303        enddo
1304
1305! Upper altitude values
1306
1307        beta(nl)=1D0/T(nl)
1308        alpha(nl)=beta(nl)*(3D0*T(nl)-4D0*T(nl-1)+T(nl-2))*Dzinv     
1309        delta(nl)=(3D0*D(nl)-4D0*D(nl-1)+D(nl-2))*Dzinv
1310
1311! Compute the gama and eps coefficients
1312! Lower altitude values
1313
1314        gama(1)=(-3D0*beta(1)+4D0*beta(2)-beta(3))*Dzinv
1315        eps(1)=(-3D0*alpha(1)+4D0*alpha(2)-alpha(3))*Dzinv
1316
1317        gama(nl-1)=(beta(nl)-beta(nl-2))*Dzinv
1318        eps(nl-1)=(alpha(nl)-alpha(nl-2))*Dzinv
1319
1320!       do l=2,nl-1
1321!       gama(l)=(beta(l+1)-beta(l-1))*Dzinv
1322!       eps(l)=(alpha(l+1)-alpha(l-1))*Dzinv
1323!       end do
1324
1325        gama(nl)=(3D0*beta(nl)-4D0*beta(nl-1)+beta(nl-2))*Dzinv
1326        eps(nl)=(3D0*alpha(nl)-4D0*alpha(nl-1)+alpha(nl-2))*Dzinv
1327
1328!       do l=1,nl
1329!       print*,'test diffparam',alpha(l),beta(l),delta(l),gama(l),eps(l)
1330!       enddo
1331!       stop
1332
1333        END
1334
1335
1336        SUBROUTINE MATCOEFF(alpha,beta,gama,delta,eps,Dad,rhoad,        &
1337     &  FacEsc,dz,dt,A,B,C,D,nl)
1338        IMPLICIT NONE
1339        INTEGER :: nl,l
1340        REAL*8, DIMENSION(nl) :: alpha,beta,gama,delta,eps,Dad,RHoad
1341        REAL*8 :: dz,dt,del1,del2,del3,FacEsc
1342        REAL*8, DIMENSION(nl) :: A,B,C,D
1343        del1=dt/2d0/dz
1344        del2=dt/dz/dz
1345        del3=dt
1346
1347! lower boundary coefficients no change
1348        A(1)=0d0
1349        B(1)=1d0
1350        C(1)=0d0
1351        D(1)=rhoAd(1)
1352
1353        do l=2,nl-1
1354        A(l)=(delta(l)+(alpha(l)+beta(l))*Dad(l))*del1-Dad(l)*del2
1355        B(l)=-(delta(l)*(alpha(l)+beta(l))+Dad(l)*(gama(l)+eps(l)))*del3
1356        B(l)=B(l)+1d0+2d0*Dad(l)*del2
1357        C(l)=-(delta(l)+(alpha(l)+beta(l))*Dad(l))*del1-Dad(l)*del2
1358        D(l)=rhoAD(l)
1359        enddo
1360
1361
1362! Upper boundary coefficients Diffusion profile
1363        C(nl)=0d0
1364        B(nl)=-1d0
1365        A(nl)=exp(-dZ*(alpha(nl)+beta(nl)))*FacEsc
1366        D(nl)=0D0
1367
1368
1369        END
1370
1371        SUBROUTINE Checkmass(X,Y,nl,nn)
1372        IMPLICIT NONE
1373
1374        INTEGER :: nl,nn
1375        REAL*8,DIMENSION(nl) :: X,Y
1376        REAL*8 Xtot,Ytot
1377
1378        Xtot=sum(X)
1379        Ytot=sum(Y)
1380
1381        if (abs((Xtot-Ytot)/Xtot) .gt. 1d-3) then
1382        print*,'no conservation for mass',Xtot,Ytot,nn
1383        endif
1384        END
1385
1386        SUBROUTINE Checkmass2(qold,qnew,P,il,nl,nn,nq)
1387        IMPLICIT NONE
1388        INTEGER :: nl,nn,l,nq,il
1389        REAL,DIMENSION(nl+1) :: P
1390        REAL*8,DIMENSION(nl,nq) :: qold,qnew
1391        REAL*8 :: DM,Mold,Mnew,g
1392        g=3.72d0
1393        DM=0d0
1394        Mold=0d0
1395        Mnew=0d0
1396        DO l=il,nl
1397        DM=DM+(qnew(l,nn)-qold(l,nn))*(dble(P(l))-dble(P(l+1)))/g
1398        Mold=Mold+qold(l,nn)*(dble(P(l))-dble(P(l+1)))/g
1399        Mnew=Mnew+qnew(l,nn)*(dble(P(l))-dble(P(l+1)))/g
1400!       print*,'l',l,qold(l,nn),qnew(l,nn),Mold,Mnew,DM,P(l),P(l+1)
1401        ENDDO
1402        IF (abs(DM/Mold) .gt. 1d-2) THEN
1403        Print*,'We dont conserve mas',nn,DM,Mold,Mnew,DM/Mold
1404        ENDIF
1405
1406        END
1407
1408        SUBROUTINE GCMGRID_P(Z,P,Q,T,Nk,Rk,qq,qnew,tt,tnew,             &
1409     &    pp,M,gc,nl,nq,nlx,ig)
1410        use tracer_mod, only: nqmx
1411        IMPLICIT NONE
1412        INTEGER :: nl,nq,nlx,il,nn,iP,ig,compteur
1413        INTEGER,DIMENSION(1) :: indP
1414        INTEGER,DIMENSION(nq) :: gc
1415        REAL*8,DIMENSION(nl) :: Z,P,T
1416        REAL*8,DIMENSION(nl,nq) :: Q,Nk,Rk
1417        REAL,DIMENSION(nqmx) :: M
1418        REAL*8,DIMENSION(nq) :: nNew
1419        REAL*8,DIMENSION(nlx) :: pp,tt,tnew
1420        REAL*8,DIMENSION(nlx) :: rhonew
1421        REAL*8,DIMENSION(nlx,nq) :: qq,qnew,rhoknew
1422        REAL*8 :: kbolt,masseU,Konst,g,Dz,facP,Hi
1423        REAL*8 :: Znew,Znew2,Pnew,Pnew2
1424        masseU=1.660538782d-27
1425        kbolt=1.3806504d-23
1426        Konst=Kbolt/masseU
1427        g=3.72d0
1428        Dz=Z(2)-Z(1)
1429        Znew=Z(nl)
1430        Znew2=Znew+dz
1431!       print*,'dz',Znew,Znew2,dz
1432        nNew(1:nq)=Nk(nl,1:nq)
1433        Pnew=P(nl)
1434
1435        do il=1,nlx
1436!       print*,'il',il,pp(il),P(1),P(nl)
1437        if (pp(il) .ge. P(1)) then
1438        qnew(il,:)=qq(il,:)
1439        tnew(il)=tt(il)
1440        endif
1441        if (pp(il) .lt. P(1)) then
1442        if (pp(il) .gt. P(nl)) then
1443        indP=maxloc(P,mask=P <  pp(il))
1444        iP=indP(1)-1
1445        if (iP .lt. 1 .or. iP .gt. nl) then
1446        print*,'danger 2',iP,nl,pp(il)
1447        endif
1448        facP=(pp(il)-P(ip))/(P(ip+1)-P(ip))
1449!       print*,'P',P(ip),P(ip+1),facP,indP,iP
1450
1451!       do nn=1,nq
1452!       qnew(il,nn)=Q(iP,nn)+
1453!     &  (Q(ip+1,nn)-Q(ip,nn))*facP
1454!       enddo
1455
1456        do nn=1,nq
1457        rhoknew(il,nn)=Rk(iP,nn)+                                       &
1458     &  (Rk(ip+1,nn)-Rk(ip,nn))*facP
1459        enddo
1460        tnew(il)=T(iP)+(T(iP+1)-T(iP))*facP
1461        rhonew(il)=sum(rhoknew(il,:))
1462        do nn=1,nq
1463        qnew(il,nn)=rhoknew(il,nn)/rhonew(il)
1464        enddo
1465
1466        else ! pp < P(nl) need to extrapolate density of each specie
1467        Pnew2=Pnew
1468
1469        compteur=0
1470        do while (pnew2 .ge. pp(il))   
1471        compteur=compteur+1
1472        do nn=1,nq
1473        Hi=Konst*T(nl)/dble(M(gc(nn)))/g
1474        Nnew(nn)=Nnew(nn)*exp(-dZ/Hi)
1475        enddo
1476        Pnew=Pnew2
1477        Pnew2=kbolt*T(nl)*sum(Nnew(:))
1478        Znew=Znew2
1479        Znew2=Znew2+Dz
1480        if (compteur .ge. 100000) then
1481        print*,'error moldiff_red infinite loop'
1482        print*,ig,il,pp(il),tt(nl),pnew2,qnew(il,:),Znew2
1483        stop
1484        endif
1485!       print*,'test',Pnew2,Znew2,Nnew(nq),pp(il)
1486        enddo
1487       
1488        facP=(pp(il)-Pnew)/(Pnew2-Pnew)
1489
1490!       do nn=1,nq
1491!       qnew(il,nn)=dble(M(gc(nn)))*Nnew(nn)
1492!     &  /sum(dble(M(gc(:)))*Nnew(:))
1493!       enddo
1494
1495        do nn=1,nq
1496        rhoknew(il,nn)=dble(M(gc(nn)))*Nnew(nn)
1497        enddo
1498        rhonew(il)=sum(rhoknew(il,:))
1499        do nn=1,nq
1500        qnew(il,nn)=rhoknew(il,nn)/rhonew(il)
1501        enddo
1502        tnew(il)=T(nl)
1503        endif
1504        endif
1505        enddo
1506
1507        END
1508
1509        SUBROUTINE GCMGRID_P2(Z,P,Q,T,Nk,Rk,qq,qnew,tt,tnew             &
1510     &   ,pp,M,gc,nl,nq,nlx,facM,ig)
1511        use tracer_mod, only: nqmx
1512        IMPLICIT NONE
1513        INTEGER :: nl,nq,nlx,il,nn,iP,ig,compteur
1514        INTEGER,DIMENSION(1) :: indP
1515        INTEGER,DIMENSION(nq) :: gc
1516        REAL*8,DIMENSION(nl) :: Z,P,T
1517        REAL*8,DIMENSION(nl,nq) :: Q,Nk,Rk
1518        REAL,DIMENSION(nqmx) :: M
1519        REAL*8,DIMENSION(nq) :: nNew
1520        REAL*8,DIMENSION(nlx) :: pp,rhonew,tt,tnew
1521        REAL*8,DIMENSION(nlx,nq) :: qq,qnew,facM,rhoknew
1522        REAL*8 :: kbolt,masseU,Konst,g,Dz,facP,Hi
1523        REAL*8 :: Znew,Znew2,Pnew,Pnew2
1524        masseU=1.660538782d-27
1525        kbolt=1.3806504d-23
1526        Konst=Kbolt/masseU
1527        g=3.72d0
1528        Dz=Z(2)-Z(1)
1529        Znew=Z(nl)
1530        Znew2=Znew+dz
1531!       print*,'dz',Znew,Znew2,dz
1532        nNew(1:nq)=Nk(nl,1:nq)
1533        Pnew=P(nl)
1534
1535        do il=1,nlx
1536!       print*,'il',il,pp(il),P(1),P(nl)
1537        if (pp(il) .ge. P(1)) then
1538        qnew(il,:)=qq(il,:)
1539        tnew(il)=tt(il)
1540        endif
1541        if (pp(il) .lt. P(1)) then
1542        if (pp(il) .gt. P(nl)) then
1543        indP=maxloc(P,mask=P <  pp(il))
1544        iP=indP(1)-1
1545        if (iP .lt. 1 .or. iP .gt. nl) then
1546        print*,'danger 3',iP,nl,pp(il)
1547        endif
1548        facP=(pp(il)-P(ip))/(P(ip+1)-P(ip))
1549!       print*,'P',P(ip),P(ip+1),facP,indP,iP
1550
1551!        do nn=1,nq
1552!        qnew(il,nn)=Q(iP,nn)+
1553!     &  (Q(ip+1,nn)-Q(ip,nn))*facP
1554!        enddo
1555
1556        do nn=1,nq
1557        rhoknew(il,nn)=(RK(iP,nn)+                                       &
1558     &  (RK(iP+1,nn)-Rk(iP,nn))*facP)*facM(il,nn)
1559        enddo
1560        tnew(il)=T(iP)+(T(ip+1)-T(iP))*facP
1561        rhonew(il)=sum(rhoknew(il,:))
1562        do nn=1,nq
1563        qnew(il,nn)=rhoknew(il,nn)/rhonew(il)
1564        enddo
1565
1566        else ! pp < P(nl) need to extrapolate density of each specie
1567        Pnew2=Pnew
1568
1569        compteur=0
1570        do while (pnew2 .ge. pp(il))   
1571        compteur=compteur+1
1572        do nn=1,nq
1573        Hi=Konst*T(nl)/dble(M(gc(nn)))/g
1574        Nnew(nn)=Nnew(nn)*exp(-dZ/Hi)
1575        enddo
1576        Pnew=Pnew2
1577        Pnew2=kbolt*T(nl)*sum(Nnew(:))
1578        Znew=Znew2
1579        Znew2=Znew2+Dz
1580        if (compteur .ge. 100000) then
1581        print*,'pb moldiff_red infinite loop'
1582        print*,ig,nl,T(nl),pnew2,qnew(il,:),Znew2
1583        stop
1584        endif
1585
1586!       print*,'test',Pnew2,Znew2,Nnew(nq),pp(il)
1587        enddo
1588       
1589        facP=(pp(il)-Pnew)/(Pnew2-Pnew)
1590
1591!        do nn=1,nq
1592!        qnew(il,nn)=dble(M(gc(nn)))*Nnew(nn)
1593!     &  /sum(dble(M(gc(:)))*Nnew(:))
1594!        enddo
1595
1596        do nn=1,nq
1597        rhoknew(il,nn)=dble(M(gc(nn)))*Nnew(nn)*FacM(il,nn)
1598        enddo
1599        rhonew(il)=sum(rhoknew(il,:))
1600        do nn=1,nq
1601        qnew(il,nn)=rhoknew(il,nn)/rhonew(il)
1602        enddo
1603        tnew(il)=T(nl)
1604
1605        endif
1606        endif
1607        enddo
1608
1609        END
1610       
Note: See TracBrowser for help on using the repository browser.