source: trunk/mesoscale/LMDZ.MARS/libf_gcm/aeronomars/moldiff.F @ 57

Last change on this file since 57 was 57, checked in by aslmd, 14 years ago

mineur LMD_MM_MARS: ajout du GCM ancienne physique, systeme maintenant complet sur SVN (ne manque que la base de donnees d'etats initiaux)

File size: 15.3 KB
Line 
1      subroutine moldiff(pplay,pplev,pt,pdt,pq,pdq,ptimestep,
2     &                   zzlay,pdteuv,pdtconduc,pdqdiff)
3
4
5      implicit none
6
7#include "dimensions.h"
8#include "dimphys.h"
9#include "comcstfi.h"
10#include "callkeys.h"
11#include "comdiurn.h"
12#include "chimiedata.h"
13#include "tracer.h"
14#include "conc.h"
15
16
17c
18c Input/Output
19c
20      real ptimestep
21      real pplay(ngridmx,nlayermx)
22      real zzlay(ngridmx,nlayermx)
23      real pplev(ngridmx,nlayermx+1)
24      real pq(ngridmx,nlayermx,nqmx)
25      real pdq(ngridmx,nlayermx,nqmx)
26      real pt(ngridmx,nlayermx)
27      real pdt(ngridmx,nlayermx)
28      real pdteuv(ngridmx,nlayermx)
29      real pdtconduc(ngridmx,nlayermx)
30      real pdqdiff(ngridmx,nlayermx,nqmx)
31c
32c Local
33c
34       real hco2(ncomptot),ho
35
36      integer ig,nz,l,n,nn
37      real del1,del2, tmean ,dalfinvdz, d
38      real hh,dcoef,dcoef1,ptfac, ntot, dens, dens2, dens3
39      real hp(nlayermx)
40      real tt(nlayermx)
41      real qq(nlayermx,ncomptot)
42      real dmmeandz(nlayermx)
43      real qnew(nlayermx,ncomptot)
44      real zlocal(nlayermx)
45      real alf(ncomptot-1,ncomptot-1)
46      real alfinv(nlayermx,ncomptot-1,ncomptot-1)
47      real indx(ncomptot-1)
48      real b(nlayermx,ncomptot-1)
49      real y(ncomptot-1,ncomptot-1)
50      real aa(nlayermx,ncomptot-1,ncomptot-1)
51      real bb(nlayermx,ncomptot-1,ncomptot-1)
52      real cc(nlayermx,ncomptot-1,ncomptot-1)
53      real atri(nlayermx-2)
54      real btri(nlayermx-2)
55      real ctri(nlayermx-2)
56      real rtri(nlayermx-2)
57      real qtri(nlayermx-2)
58      real alfdiag(ncomptot-1)
59      real wi(ncomptot), flux(ncomptot), pote
60
61      integer i_co2, i_co, i_o2, i_h2, i_h2o, i_h2o2,
62     $        i_o1d, i_o, i_h, i_oh, i_ho2, i_n2, i_o3, i_ar
63      integer g_co2, g_co, g_o2, g_h2, g_h2o, g_h2o2,
64     $        g_o1d, g_o, g_h, g_oh, g_ho2, g_o3, g_n2, g_ar
65      integer gcmind(ncomptot)
66      integer ierr
67
68      logical firstcall
69      real abfac(ncomptot)
70      real dij(ncomptot,ncomptot)
71      save firstcall
72      save dij
73      data firstcall /.true./
74
75      if (firstcall) then
76        call moldiffcoeff(dij)
77        print*,'MOLDIFF  EXO'
78        firstcall= .false.
79      endif
80
81cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
82c     tracer numbering in the gcm
83cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
84c
85      g_co2      =  nqchem_min
86      g_co       =  nqchem_min + 1
87      g_o        =  nqchem_min + 2
88      g_o1d      =  nqchem_min + 3
89      g_o2       =  nqchem_min + 4
90      g_o3       =  nqchem_min + 5
91      g_h        =  nqchem_min + 6
92      g_h2       =  nqchem_min + 7
93      g_oh       =  nqchem_min + 8
94      g_ho2      =  nqchem_min + 9
95      g_h2o2     =  nqchem_min + 10
96      g_n2       =  nqchem_min + 11
97      g_ar       =  nqchem_min + 12
98      g_h2o      =  nqmx
99
100cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
101c     tracer numbering in the molecular diffusion
102cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
103c  Atomic oxygen must always be the LAST species of the list as
104c it is the dominant species at high altitudes. 
105
106      i_co   = 1
107      i_n2   = 2
108      i_o2   = 3
109      i_co2  = 4
110      i_h2   = 5
111      i_h    = 6
112      i_oh   = 7
113      i_ho2  = 8
114      i_h2o  = 9
115      i_h2o2 = 10
116      i_o1d  = 11
117      i_o3   = 12
118      i_ar   = 13
119      i_o    = 14
120
121cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
122c     array to relate local indexes to gcm indexes
123cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
124
125      gcmind(i_co)      =       g_co
126      gcmind(i_n2)      =       g_n2
127      gcmind(i_o2)      =       g_o2
128      gcmind(i_co2)     =       g_co2
129      gcmind(i_h2)      =       g_h2
130      gcmind(i_h)       =       g_h
131      gcmind(i_oh)      =       g_oh
132      gcmind(i_ho2)     =       g_ho2
133      gcmind(i_h2o)     =       g_h2o
134      gcmind(i_h2o2)=   g_h2o2
135      gcmind(i_o1d) =   g_o1d
136      gcmind(i_o3)  =   g_o3
137      gcmind(i_o)   =   g_o
138      gcmind(i_ar)   =  g_ar
139c
140cccccccccccccccccccccccccccccccccccccccccccccccccccccccc
141
142      nz=nlayermx
143
144      do ig=1,ngridmx
145
146        do l=2,nz-1
147          tt(l)=pt(ig,l)+pdt(ig,l)*ptimestep
148     &                  +pdteuv(ig,l)*ptimestep
149     &                  +pdtconduc(ig,l)*ptimestep
150
151          do nn=1,ncomptot
152            qq(l,nn)=pq(ig,l,gcmind(nn))+pdq(ig,l,gcmind(nn))*ptimestep
153            qq(l,nn)=max(qq(l,nn),1.e-30)
154          enddo
155          hp(l)=-log(pplay(ig,l+1)/pplay(ig,l-1))
156          dmmeandz(l)=(mmean(ig,l+1)-mmean(ig,l-1))/hp(l)
157        enddo
158
159        tt(1)=pt(ig,1)  +pdt(ig,1)*ptimestep
160     &                  +pdteuv(ig,1)*ptimestep
161     &                  +pdtconduc(ig,1)*ptimestep
162        tt(nz)=pt(ig,nz)+pdt(ig,nz)*ptimestep
163     &                  +pdteuv(ig,nz)*ptimestep
164     &                  +pdtconduc(ig,nz)*ptimestep
165
166        do nn=1,ncomptot
167          qq(1,nn)=pq(ig,1,gcmind(nn))+pdq(ig,1,gcmind(nn))*ptimestep
168          qq(nz,nn)=pq(ig,nz,gcmind(nn))+pdq(ig,nz,gcmind(nn))*ptimestep
169          qq(1,nn)=max(qq(1,nn),1.e-30)
170          qq(nz,nn)=max(qq(nz,nn),1.e-30)
171        enddo
172        hp(1)=-log(pplay(ig,2)/pplay(ig,1))
173        dmmeandz(1)=(-3.*mmean(ig,1)+4.*mmean(ig,2)
174     &               -mmean(ig,3))/(2.*hp(1))
175        hp(nz)=-log(pplay(ig,nz)/pplay(ig,nz-1))
176        dmmeandz(nz)=(3.*mmean(ig,nz)-4.*mmean(ig,nz-1)
177     &                +mmean(ig,nz-2))/(2.*hp(nz))
178c
179c Setting-up matrix of alfa coefficients from Dickinson and Ridley 1972
180c
181      do l=1,nz
182       if(abs(dmmeandz(l)) .lt. 1.e-5) dmmeandz(l)=0.0
183        hh=rnew(ig,l)*tt(l)/g
184        ptfac=(1.e5/pplay(ig,l))*(tt(l)/273)**1.75
185        ntot=pplay(ig,l)/(1.381e-23*tt(l))                      ! in #/m3
186
187        do nn=1,ncomptot-1            ! rows
188          alfdiag(nn)=0.
189          dcoef1=dij(nn,i_o)*ptfac
190          do n=1,ncomptot-1           ! columns
191            y(nn,n)=0.
192            dcoef=dij(nn,n)*ptfac
193            alf(nn,n)=qq(l,nn)/ntot/1.66e-27
194     &         *(1./(mmol(gcmind(n))*dcoef)-1./(mmol(g_o)*dcoef1))
195            alfdiag(nn)=alfdiag(nn)
196     &       +(1./(mmol(gcmind(n))*dcoef)-1./(mmol(g_o)*dcoef1))
197     &        *qq(l,n)
198          enddo
199          dcoef=dij(nn,nn)*ptfac
200          alfdiag(nn)=alfdiag(nn)
201     &       -(1./(mmol(gcmind(nn))*dcoef)-1./(mmol(g_o)*dcoef1))
202     &          *qq(l,nn)
203          alf(nn,nn)=-(alfdiag(nn)
204     &                 +1./(mmol(g_o)*dcoef1))/ntot/1.66e-27
205          y(nn,nn)=1.
206          b(l,nn)=-(dmmeandz(l)/mmean(ig,l)
207     &              +mmol(gcmind(nn))/mmean(ig,l)-1.)
208        enddo
209
210c
211c Inverting the alfa matrix
212c
213        call ludcmp(alf,ncomptot-1,ncomptot-1,indx,d,ierr)
214
215c       TEMPORAIRE *****************************
216        if (ierr.ne.0) then
217            write(*,*) 'In moldiff: Problem in LUDCMP with matrix alf'
218            write(*,*) 'Singular matrix ?'
219c           write(*,*) 'Matrix alf = ', alf
220            write(*,*) 'ig, l=',ig, l
221            write(*,*) 'No molecular diffusion this time !'
222            call zerophys(ngridmx*nlayermx*nqmx,pdqdiff)
223            return
224c           stop
225        end if
226c       *******************************************
227        do n=1,ncomptot-1
228          call lubksb(alf,ncomptot-1,ncomptot-1,indx,y(1,n))
229          do nn=1,ncomptot-1
230            alfinv(l,nn,n)=y(nn,n)/hh
231          enddo
232        enddo
233      enddo
234
235c
236c Calculating coefficients of the system
237c
238
239c      zlocal(1)=-log(pplay(ig,1)/pplev(ig,1))* Rnew(ig,1)*tt(1)/g
240      zlocal(1)=zzlay(ig,1)
241
242      do l=2,nz-1
243        del1=hp(l)*pplay(ig,l)/(g*ptimestep)
244        del2=(hp(l)/2)**2*pplay(ig,l)/(g*ptimestep)
245        do nn=1,ncomptot-1
246          do n=1,ncomptot-1
247            dalfinvdz=(alfinv(l+1,nn,n)-alfinv(l-1,nn,n))/hp(l)
248            aa(l,nn,n)=-dalfinvdz/del1+alfinv(l,nn,n)/del2
249     &                +alfinv(l-1,nn,n)*b(l-1,n)/del1   
250            bb(l,nn,n)=-2.*alfinv(l,nn,n)/del2
251            cc(l,nn,n)=dalfinvdz/del1+alfinv(l,nn,n)/del2
252     &                -alfinv(l+1,nn,n)*b(l+1,n)/del1   
253          enddo
254        enddo
255
256c        tmean=tt(l)
257c        if(tt(l).ne.tt(l-1))
258c     &        tmean=(tt(l)-tt(l-1))/log(tt(l)/tt(l-1))
259c        zlocal(l)= zlocal(l-1)
260c     &         -log(pplay(ig,l)/pplay(ig,l-1))*rnew(ig,l)*tmean/g
261      zlocal(l)=zzlay(ig,l)
262      enddo
263
264c      zlocal(nz)= zlocal(nz-1)
265c     &         -log(pplay(ig,nz)/pplay(ig,nz-1))*rnew(ig,nz)*tmean/g
266      zlocal(nz)=zzlay(ig,nz)
267       
268ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
269c Escape velocity from Jeans equation for the flux of H and H2
270c (Hunten 1973, eq. 5)
271     
272      do n=1,ncomptot
273        wi(n)=1.
274        flux(n)=0.
275        abfac(n)=1.
276      enddo
277
278       dens=pplay(ig,nz)/(rnew(ig,nz)*tt(nz))
279c
280c For H:
281c
282       pote=(3398000.+zlocal(nz))/
283     &         (1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h)*g))
284       wi(i_h)=sqrt(2.*1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h)))
285     &             /(2.*sqrt(3.1415))*(1.+pote)*exp(-pote)
286       flux(i_h)=qq(nz,i_h)*dens/(1.6605e-27*mmol(g_h))*wi(i_h)
287       flux(i_h)=flux(i_h)*1.6606e-27
288       abfac(i_h)=0.
289c
290c For H2:
291c
292       pote=(3398000.+zlocal(nz))/
293     &         (1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h2)*g))
294       wi(i_h2)=sqrt(2.*1.381e-23*tt(nz)/(1.6605e-27*mmol(g_h2)))
295     &              /(2.*sqrt(3.1415))*(1.+pote)*exp(-pote)
296       flux(i_h2)=qq(nz,i_h2)*dens/(1.6605e-27*mmol(g_h2))*wi(i_h2)
297       flux(i_h2)=flux(i_h2)*1.6606e-27
298       abfac(i_h2)=0.
299
300c ********* TEMPORAIRE : no escape for h and h2
301c     do n=1,ncomptot
302c       wi(n)=1.
303c       flux(n)=0.
304c       abfac(n)=1.
305c     enddo
306c ********************************************
307
308
309ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
310
311c
312c Setting coefficients for tridiagonal matrix and solving the system
313c
314
315      do nn=1,ncomptot-1
316        do l=2,nz-1
317          atri(l-1)=aa(l,nn,nn)
318          btri(l-1)=bb(l,nn,nn)+1.
319          ctri(l-1)=cc(l,nn,nn)
320          rtri(l-1)=qq(l,nn)
321          do n=1,ncomptot-1
322            rtri(l-1)=rtri(l-1)-(aa(l,nn,n)*qq(l-1,n)
323     &                           +bb(l,nn,n)*qq(l,n)
324     &                           +cc(l,nn,n)*qq(l+1,n))
325          enddo
326          rtri(l-1)=rtri(l-1)+(aa(l,nn,nn)*qq(l-1,nn)
327     &                         +bb(l,nn,nn)*qq(l,nn)
328     &                         +cc(l,nn,nn)*qq(l+1,nn))
329        enddo
330
331c
332c Boundary conditions:
333c       Escape flux for H and H2 at top
334c       Diffusive equilibrium for the other species at top
335c       Perfect mixing for all at bottom
336c
337
338        rtri(nz-2)=rtri(nz-2)
339     &             -ctri(nz-2)*flux(nn)*mmol(gcmind(nn))/(dens*wi(nn))
340
341        atri(nz-2)=atri(nz-2)
342     &             -abfac(nn)*ctri(nz-2)/(3.-2.*hp(nz)*b(nz,nn))
343        btri(nz-2)=btri(nz-2)
344     &             +abfac(nn)*4.*ctri(nz-2)/(3.-2.*hp(nz)*b(nz,nn))
345
346c        rtri(1)=rtri(1)-atri(1)*qq(1,nn)
347        btri(1)=btri(1)+atri(1)
348
349        call tridag(atri,btri,ctri,rtri,qtri,nz-2)
350
351        do l=2,nz-1
352c          qnew(l,nn)=qtri(l-1)
353          qnew(l,nn)=max(qtri(l-1),1.e-30)
354        enddo
355
356        qnew(nz,nn)=flux(nn)*mmol(gcmind(nn))/(dens*wi(nn))
357     &               +abfac(nn)*(4.*qnew(nz-1,nn)-qnew(nz-2,nn))
358     &                /(3.-2.*hp(nz)*b(nz,nn))
359c        qnew(1,nn)=qq(1,nn)
360        qnew(1,nn)=qnew(2,nn)
361         
362        qnew(nz,nn)=max(qnew(nz,nn),1.e-30)
363        qnew(1,nn)=max(qnew(1,nn),1.e-30)
364
365      enddo     ! loop on species
366
367      DO l=1,nz
368        if(zlocal(l).gt.65000.)then
369        pdqdiff(ig,l,g_o)=0.
370        do n=1,ncomptot-1
371          pdqdiff(ig,l,gcmind(n))=(qnew(l,n)-qq(l,n))
372          pdqdiff(ig,l,g_o)=pdqdiff(ig,l,g_o)-(qnew(l,n)-qq(l,n))
373          pdqdiff(ig,l,gcmind(n))=pdqdiff(ig,l,gcmind(n))/ptimestep
374        enddo
375                pdqdiff(ig,l,g_o)=pdqdiff(ig,l,g_o)/ptimestep
376        endif
377      ENDDO
378
379c      do l=2,nz
380c        do n=1,ncomptot-1
381c          hco2(n)=qnew(l,n)*pplay(ig,l)/(rnew(ig,l)*tt(l)) /
382c     &      (qnew(l-1,n)*pplay(ig,l-1)/(rnew(ig,l-1)*tt(l-1)))
383c          hco2(n)=-(zlocal(l)-zlocal(l-1))/log(hco2(n))/1000.
384c        enddo
385c        write(225,*),l,pt(1,l),(hco2(n),n=1,6)
386c        write(226,*),l,pt(1,l),(hco2(n),n=7,12)
387c      enddo
388
389      enddo             ! ig loop
390
391      return
392      end
393
394c    ********************************************************************
395c    ********************************************************************
396c    ********************************************************************
397 
398      subroutine tridag(a,b,c,r,u,n)
399      parameter (nmax=100)
400c      dimension gam(nmax),a(n),b(n),c(n),r(n),u(n)
401      real gam(nmax),a(n),b(n),c(n),r(n),u(n)
402      if(b(1).eq.0.)pause
403      bet=b(1)
404      u(1)=r(1)/bet
405      do 11 j=2,n
406        gam(j)=c(j-1)/bet
407        bet=b(j)-a(j)*gam(j)
408        if(bet.eq.0.)pause
409        u(j)=(r(j)-a(j)*u(j-1))/bet
41011    continue
411      do 12 j=n-1,1,-1
412        u(j)=u(j)-gam(j+1)*u(j+1)
41312    continue
414      return
415      end
416
417c    ********************************************************************
418c    ********************************************************************
419c    ********************************************************************
420
421      SUBROUTINE LUBKSB(A,N,NP,INDX,B)
422
423      implicit none
424
425      integer i,j,n,np,ii,ll
426      real sum
427      real a(np,np),indx(np),b(np)
428
429c      DIMENSION A(NP,NP),INDX(N),B(N)
430      II=0
431      DO 12 I=1,N
432        LL=INDX(I)
433        SUM=B(LL)
434        B(LL)=B(I)
435        IF (II.NE.0)THEN
436          DO 11 J=II,I-1
437            SUM=SUM-A(I,J)*B(J)
43811        CONTINUE
439        ELSE IF (SUM.NE.0.) THEN
440          II=I
441        ENDIF
442        B(I)=SUM
44312    CONTINUE
444      DO 14 I=N,1,-1
445        SUM=B(I)
446        IF(I.LT.N)THEN
447          DO 13 J=I+1,N
448            SUM=SUM-A(I,J)*B(J)
44913        CONTINUE
450        ENDIF
451        B(I)=SUM/A(I,I)
45214    CONTINUE
453      RETURN
454      END
455
456c    ********************************************************************
457c    ********************************************************************
458c    ********************************************************************
459
460      SUBROUTINE LUDCMP(A,N,NP,INDX,D,ierr)
461
462      implicit none
463
464      integer n,np,nmax,i,j,k,imax
465      real d,tiny,aamax
466      real a(np,np),indx(np)
467      integer ierr  ! error =0 if OK, =1 if problem
468
469      PARAMETER (NMAX=100,TINY=1.0E-20)
470c      DIMENSION A(NP,NP),INDX(N),VV(NMAX)
471      real sum,vv(nmax),dum
472
473      D=1.
474      DO 12 I=1,N
475        AAMAX=0.
476        DO 11 J=1,N
477          IF (ABS(A(I,J)).GT.AAMAX) AAMAX=ABS(A(I,J))
47811      CONTINUE
479        IF (AAMAX.EQ.0.) then
480            write(*,*) 'In moldiff: Problem in LUDCMP with matrix A'
481            write(*,*) 'Singular matrix ?'
482c           write(*,*) 'Matrix A = ', A
483c           TO DEBUG :
484            ierr =1
485            return
486c           stop
487        END IF
488
489        VV(I)=1./AAMAX
49012    CONTINUE
491      DO 19 J=1,N
492        IF (J.GT.1) THEN
493          DO 14 I=1,J-1
494            SUM=A(I,J)
495            IF (I.GT.1)THEN
496              DO 13 K=1,I-1
497                SUM=SUM-A(I,K)*A(K,J)
49813            CONTINUE
499              A(I,J)=SUM
500            ENDIF
50114        CONTINUE
502        ENDIF
503        AAMAX=0.
504        DO 16 I=J,N
505          SUM=A(I,J)
506          IF (J.GT.1)THEN
507            DO 15 K=1,J-1
508              SUM=SUM-A(I,K)*A(K,J)
50915          CONTINUE
510            A(I,J)=SUM
511          ENDIF
512          DUM=VV(I)*ABS(SUM)
513          IF (DUM.GE.AAMAX) THEN
514            IMAX=I
515            AAMAX=DUM
516          ENDIF
51716      CONTINUE
518        IF (J.NE.IMAX)THEN
519          DO 17 K=1,N
520            DUM=A(IMAX,K)
521            A(IMAX,K)=A(J,K)
522            A(J,K)=DUM
52317        CONTINUE
524          D=-D
525          VV(IMAX)=VV(J)
526        ENDIF
527        INDX(J)=IMAX
528        IF(J.NE.N)THEN
529          IF(A(J,J).EQ.0.)A(J,J)=TINY
530          DUM=1./A(J,J)
531          DO 18 I=J+1,N
532            A(I,J)=A(I,J)*DUM
53318        CONTINUE
534        ENDIF
53519    CONTINUE
536      IF(A(N,N).EQ.0.)A(N,N)=TINY
537      ierr =0
538      RETURN
539      END
540
Note: See TracBrowser for help on using the repository browser.