source: LMDZ4/trunk/libf/phylmd/yamada4.F @ 3839

Last change on this file since 3839 was 1403, checked in by Laurent Fairhead, 14 years ago

Merged LMDZ4V5.0-dev branch changes r1292:r1399 to trunk.

Validation:
Validation consisted in compiling the HEAD revision of the trunk,
LMDZ4V5.0-dev branch and the merged sources and running different
configurations on local and SX8 machines comparing results.

Local machine: bench configuration, 32x24x11, gfortran

  • IPSLCM5A configuration (comparison between trunk and merged sources):
    • numerical convergence on dynamical fields over 3 days
    • start files are equivalent (except for RN and PB fields)
    • daily history files equivalent
  • MH07 configuration, new physics package (comparison between LMDZ4V5.0-dev branch and merged sources):
    • numerical convergence on dynamical fields over 3 days
    • start files are equivalent (except for RN and PB fields)
    • daily history files equivalent

SX8 machine (brodie), 96x95x39 on 4 processors:

  • IPSLCM5A configuration:
    • start files are equivalent (except for RN and PB fields)
    • monthly history files equivalent
  • MH07 configuration:
    • start files are equivalent (except for RN and PB fields)
    • monthly history files equivalent

Changes to the makegcm and create_make_gcm scripts to take into account
main programs in F90 files


Fusion de la branche LMDZ4V5.0-dev (r1292:r1399) au tronc principal

Validation:
La validation a consisté à compiler la HEAD de le trunk et de la banche
LMDZ4V5.0-dev et les sources fusionnées et de faire tourner le modéle selon
différentes configurations en local et sur SX8 et de comparer les résultats

En local: 32x24x11, config bench/gfortran

  • pour une config IPSLCM5A (comparaison tronc/fusion):
    • convergence numérique sur les champs dynamiques après 3 jours
    • restart et restartphy égaux (à part sur RN et Pb)
    • fichiers histoire égaux
  • pour une config nlle physique (MH07) (comparaison LMDZ4v5.0-dev/fusion):
    • convergence numérique sur les champs dynamiques après 3 jours
    • restart et restartphy égaux
    • fichiers histoire équivalents

Sur brodie, 96x95x39 sur 4 proc:

  • pour une config IPSLCM5A:
    • restart et restartphy égaux (à part sur RN et PB)
    • pas de différence dans les fichiers histmth.nc
  • pour une config MH07
    • restart et restartphy égaux (à part sur RN et PB)
    • pas de différence dans les fichiers histmth.nc

Changement sur makegcm et create_make-gcm pour pouvoir prendre en compte des
programmes principaux en *F90

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 19.9 KB
RevLine 
[541]1!
2! $Header$
3!
4      SUBROUTINE yamada4(ngrid,dt,g,rconst,plev,temp
5     s   ,zlev,zlay,u,v,teta,cd,q2,km,kn,kq,ustar
6     s   ,iflag_pbl)
[766]7      use dimphy
[541]8      IMPLICIT NONE
[938]9#include "iniprint.h"
[541]10c.......................................................................
[766]11cym#include "dimensions.h"
12cym#include "dimphy.h"
[541]13c.......................................................................
14c
15c dt : pas de temps
16c g  : g
17c zlev : altitude a chaque niveau (interface inferieure de la couche
18c        de meme indice)
19c zlay : altitude au centre de chaque couche
20c u,v : vitesse au centre de chaque couche
21c       (en entree : la valeur au debut du pas de temps)
22c teta : temperature potentielle au centre de chaque couche
23c        (en entree : la valeur au debut du pas de temps)
24c cd : cdrag
25c      (en entree : la valeur au debut du pas de temps)
26c q2 : $q^2$ au bas de chaque couche
27c      (en entree : la valeur au debut du pas de temps)
28c      (en sortie : la valeur a la fin du pas de temps)
29c km : diffusivite turbulente de quantite de mouvement (au bas de chaque
30c      couche)
31c      (en sortie : la valeur a la fin du pas de temps)
32c kn : diffusivite turbulente des scalaires (au bas de chaque couche)
33c      (en sortie : la valeur a la fin du pas de temps)
34c
35c  iflag_pbl doit valoir entre 6 et 9
36c      l=6, on prend  systematiquement une longueur d'equilibre
37c    iflag_pbl=6 : MY 2.0
38c    iflag_pbl=7 : MY 2.0.Fournier
39c    iflag_pbl=8 : MY 2.5
[1403]40c    iflag_pbl>=9 : MY 2.5 avec diffusion verticale
[541]41
42c.......................................................................
43      REAL dt,g,rconst
44      real plev(klon,klev+1),temp(klon,klev)
45      real ustar(klon)
46      real kmin,qmin,pblhmin(klon),coriol(klon)
47      REAL zlev(klon,klev+1)
48      REAL zlay(klon,klev)
49      REAL u(klon,klev)
50      REAL v(klon,klev)
51      REAL teta(klon,klev)
52      REAL cd(klon)
53      REAL q2(klon,klev+1),qpre
54      REAL unsdz(klon,klev)
55      REAL unsdzdec(klon,klev+1)
56
57      REAL km(klon,klev+1)
58      REAL kmpre(klon,klev+1),tmp2
59      REAL mpre(klon,klev+1)
60      REAL kn(klon,klev+1)
61      REAL kq(klon,klev+1)
62      real ff(klon,klev+1),delta(klon,klev+1)
63      real aa(klon,klev+1),aa0,aa1
64      integer iflag_pbl,ngrid
65
66
67      integer nlay,nlev
68
69      logical first
70      integer ipas
71      save first,ipas
[926]72cFH/IM     data first,ipas/.true.,0/
73      data first,ipas/.false.,0/
[766]74c$OMP THREADPRIVATE( first,ipas)
[541]75
76      integer ig,k
77
78
79      real ri,zrif,zalpha,zsm,zsn
80      real rif(klon,klev+1),sm(klon,klev+1),alpha(klon,klev)
81
82      real m2(klon,klev+1),dz(klon,klev+1),zq,n2(klon,klev+1)
83      real dtetadz(klon,klev+1)
84      real m2cstat,mcstat,kmcstat
[766]85      real l(klon,klev+1)
86      real,allocatable,save :: l0(:)
87c$OMP THREADPRIVATE(l0)     
[541]88      real sq(klon),sqz(klon),zz(klon,klev+1)
89      integer iter
90
91      real ric,rifc,b1,kap
92      save ric,rifc,b1,kap
93      data ric,rifc,b1,kap/0.195,0.191,16.6,0.4/
[766]94c$OMP THREADPRIVATE(ric,rifc,b1,kap)
[541]95      real frif,falpha,fsm
96      real fl,zzz,zl0,zq2,zn2
97
[1403]98      real rino(klon,klev+1),smyam(klon,klev),styam(klon,klev)
99     s  ,lyam(klon,klev),knyam(klon,klev)
100     s  ,w2yam(klon,klev),t2yam(klon,klev)
[766]101      logical,save :: firstcall=.true.
102c$OMP THREADPRIVATE(firstcall)       
[541]103      frif(ri)=0.6588*(ri+0.1776-sqrt(ri*ri-0.3221*ri+0.03156))
104      falpha(ri)=1.318*(0.2231-ri)/(0.2341-ri)
105      fsm(ri)=1.96*(0.1912-ri)*(0.2341-ri)/((1.-ri)*(0.2231-ri))
106      fl(zzz,zl0,zq2,zn2)=
107     s     max(min(l0(ig)*kap*zlev(ig,k)/(kap*zlev(ig,k)+l0(ig))
108     s     ,0.5*sqrt(q2(ig,k))/sqrt(max(n2(ig,k),1.e-10))) ,1.)
109
[766]110
111      nlay=klev
112      nlev=klev+1
113     
114      if (firstcall) then
115        allocate(l0(klon))
116        firstcall=.false.
117      endif
118
119
[1403]120      if (.not.(iflag_pbl.ge.6.and.iflag_pbl.le.10)) then
[541]121           stop'probleme de coherence dans appel a MY'
122      endif
123
124      ipas=ipas+1
125      if (0.eq.1.and.first) then
126      do ig=1,1000
127         ri=(ig-800.)/500.
128         if (ri.lt.ric) then
129            zrif=frif(ri)
130         else
131            zrif=rifc
132         endif
133         if(zrif.lt.0.16) then
134            zalpha=falpha(zrif)
135            zsm=fsm(zrif)
136         else
137            zalpha=1.12
138            zsm=0.085
139         endif
140c     print*,ri,rif,zalpha,zsm
141      enddo
142      endif
143
144c.......................................................................
145c  les increments verticaux
146c.......................................................................
147c
148c!!!!! allerte !!!!!c
149c!!!!! zlev n'est pas declare a nlev !!!!!c
150c!!!!! ---->
151                                                      DO ig=1,ngrid
152            zlev(ig,nlev)=zlay(ig,nlay)
153     &             +( zlay(ig,nlay) - zlev(ig,nlev-1) )
154                                                      ENDDO
155c!!!!! <----
156c!!!!! allerte !!!!!c
157c
158      DO k=1,nlay
159                                                      DO ig=1,ngrid
160        unsdz(ig,k)=1.E+0/(zlev(ig,k+1)-zlev(ig,k))
161                                                      ENDDO
162      ENDDO
163                                                      DO ig=1,ngrid
164      unsdzdec(ig,1)=1.E+0/(zlay(ig,1)-zlev(ig,1))
165                                                      ENDDO
166      DO k=2,nlay
167                                                      DO ig=1,ngrid
168        unsdzdec(ig,k)=1.E+0/(zlay(ig,k)-zlay(ig,k-1))
169                                                     ENDDO
170      ENDDO
171                                                      DO ig=1,ngrid
172      unsdzdec(ig,nlay+1)=1.E+0/(zlev(ig,nlay+1)-zlay(ig,nlay))
173                                                     ENDDO
174c
175c.......................................................................
176
177      do k=2,klev
178                                                          do ig=1,ngrid
179         dz(ig,k)=zlay(ig,k)-zlay(ig,k-1)
180         m2(ig,k)=((u(ig,k)-u(ig,k-1))**2+(v(ig,k)-v(ig,k-1))**2)
181     s             /(dz(ig,k)*dz(ig,k))
182         dtetadz(ig,k)=(teta(ig,k)-teta(ig,k-1))/dz(ig,k)
183         n2(ig,k)=g*2.*dtetadz(ig,k)/(teta(ig,k-1)+teta(ig,k))
184c        n2(ig,k)=0.
185         ri=n2(ig,k)/max(m2(ig,k),1.e-10)
186         if (ri.lt.ric) then
187            rif(ig,k)=frif(ri)
188         else
189            rif(ig,k)=rifc
190         endif
191         if(rif(ig,k).lt.0.16) then
192            alpha(ig,k)=falpha(rif(ig,k))
193            sm(ig,k)=fsm(rif(ig,k))
194         else
195            alpha(ig,k)=1.12
196            sm(ig,k)=0.085
197         endif
198         zz(ig,k)=b1*m2(ig,k)*(1.-rif(ig,k))*sm(ig,k)
199c     print*,'RIF L=',k,rif(ig,k),ri*alpha(ig,k)
200
201
202                                                          enddo
203      enddo
204
205
206c====================================================================
207c   Au premier appel, on determine l et q2 de facon iterative.
208c iterration pour determiner la longueur de melange
209
210
211      if (first.or.iflag_pbl.eq.6) then
212                                                          do ig=1,ngrid
213      l0(ig)=10.
214                                                          enddo
215      do k=2,klev-1
216                                                          do ig=1,ngrid
217        l(ig,k)=l0(ig)*kap*zlev(ig,k)/(kap*zlev(ig,k)+l0(ig))
218                                                          enddo
219      enddo
220
221      do iter=1,10
222                                                          do ig=1,ngrid
223         sq(ig)=1.e-10
224         sqz(ig)=1.e-10
225                                                          enddo
226         do k=2,klev-1
227                                                          do ig=1,ngrid
228           q2(ig,k)=l(ig,k)**2*zz(ig,k)
229           l(ig,k)=fl(zlev(ig,k),l0(ig),q2(ig,k),n2(ig,k))
230           zq=sqrt(q2(ig,k))
231           sqz(ig)=sqz(ig)+zq*zlev(ig,k)*(zlay(ig,k)-zlay(ig,k-1))
232           sq(ig)=sq(ig)+zq*(zlay(ig,k)-zlay(ig,k-1))
233                                                          enddo
234         enddo
235                                                          do ig=1,ngrid
236         l0(ig)=0.2*sqz(ig)/sq(ig)
237c        l0(ig)=30.
238                                                          enddo
239c      print*,'ITER=',iter,'  L0=',l0
240
241      enddo
242
243c     print*,'Fin de l initialisation de q2 et l0'
244
245      endif ! first
246
247c====================================================================
248c  Calcul de la longueur de melange.
249c====================================================================
250
251c   Mise a jour de l0
252                                                          do ig=1,ngrid
253      sq(ig)=1.e-10
254      sqz(ig)=1.e-10
255                                                          enddo
256      do k=2,klev-1
257                                                          do ig=1,ngrid
258        zq=sqrt(q2(ig,k))
259        sqz(ig)=sqz(ig)+zq*zlev(ig,k)*(zlay(ig,k)-zlay(ig,k-1))
260        sq(ig)=sq(ig)+zq*(zlay(ig,k)-zlay(ig,k-1))
261                                                          enddo
262      enddo
263                                                          do ig=1,ngrid
264      l0(ig)=0.2*sqz(ig)/sq(ig)
265c        l0(ig)=30.
266                                                          enddo
267c      print*,'ITER=',iter,'  L0=',l0
268c   calcul de l(z)
269      do k=2,klev
270                                                          do ig=1,ngrid
271         l(ig,k)=fl(zlev(ig,k),l0(ig),q2(ig,k),n2(ig,k))
272         if(first) then
273           q2(ig,k)=l(ig,k)**2*zz(ig,k)
274         endif
275                                                          enddo
276      enddo
277
278c====================================================================
279c   Yamada 2.0
280c====================================================================
281      if (iflag_pbl.eq.6) then
282
283      do k=2,klev
284                                                          do ig=1,ngrid
285         q2(ig,k)=l(ig,k)**2*zz(ig,k)
286                                                          enddo
287      enddo
288
289
290      else if (iflag_pbl.eq.7) then
291c====================================================================
292c   Yamada 2.Fournier
293c====================================================================
294
295c  Calcul de l,  km, au pas precedent
296      do k=2,klev
297                                                          do ig=1,ngrid
298c        print*,'SMML=',sm(ig,k),l(ig,k)
299         delta(ig,k)=q2(ig,k)/(l(ig,k)**2*sm(ig,k))
300         kmpre(ig,k)=l(ig,k)*sqrt(q2(ig,k))*sm(ig,k)
301         mpre(ig,k)=sqrt(m2(ig,k))
302c        print*,'0L=',k,l(ig,k),delta(ig,k),km(ig,k)
303                                                          enddo
304      enddo
305
306      do k=2,klev-1
307                                                          do ig=1,ngrid
308        m2cstat=max(alpha(ig,k)*n2(ig,k)+delta(ig,k)/b1,1.e-12)
309        mcstat=sqrt(m2cstat)
310
311c        print*,'M2 L=',k,mpre(ig,k),mcstat
312c
313c  -----{puis on ecrit la valeur de q qui annule l'equation de m
314c        supposee en q3}
315c
316        IF (k.eq.2) THEN
317          kmcstat=1.E+0 / mcstat
318     &    *( unsdz(ig,k)*kmpre(ig,k+1)
319     &                        *mpre(ig,k+1)
320     &      +unsdz(ig,k-1)
321     &              *cd(ig)
322     &              *( sqrt(u(ig,3)**2+v(ig,3)**2)
323     &                -mcstat/unsdzdec(ig,k)
324     &                -mpre(ig,k+1)/unsdzdec(ig,k+1) )**2)
325     &      /( unsdz(ig,k)+unsdz(ig,k-1) )
326        ELSE
327          kmcstat=1.E+0 / mcstat
328     &    *( unsdz(ig,k)*kmpre(ig,k+1)
329     &                        *mpre(ig,k+1)
330     &      +unsdz(ig,k-1)*kmpre(ig,k-1)
331     &                          *mpre(ig,k-1) )
332     &      /( unsdz(ig,k)+unsdz(ig,k-1) )
333        ENDIF
334c       print*,'T2 L=',k,tmp2
335        tmp2=kmcstat
336     &      /( sm(ig,k)/q2(ig,k) )
337     &      /l(ig,k)
338        q2(ig,k)=max(tmp2,1.e-12)**(2./3.)
339c       print*,'Q2 L=',k,q2(ig,k)
340c
341                                                          enddo
342      enddo
343
344      else if (iflag_pbl.ge.8) then
345c====================================================================
346c   Yamada 2.5 a la Didi
347c====================================================================
348
349
350c  Calcul de l,  km, au pas precedent
351      do k=2,klev
352                                                          do ig=1,ngrid
353c        print*,'SMML=',sm(ig,k),l(ig,k)
354         delta(ig,k)=q2(ig,k)/(l(ig,k)**2*sm(ig,k))
355         if (delta(ig,k).lt.1.e-20) then
356c     print*,'ATTENTION   L=',k,'   Delta=',delta(ig,k)
357            delta(ig,k)=1.e-20
358         endif
359         km(ig,k)=l(ig,k)*sqrt(q2(ig,k))*sm(ig,k)
360         aa0=
361     s   (m2(ig,k)-alpha(ig,k)*n2(ig,k)-delta(ig,k)/b1)
362         aa1=
363     s   (m2(ig,k)*(1.-rif(ig,k))-delta(ig,k)/b1)
364c abder      print*,'AA L=',k,aa0,aa1,aa1/max(m2(ig,k),1.e-20)
365         aa(ig,k)=aa1*dt/(delta(ig,k)*l(ig,k))
366c     print*,'0L=',k,l(ig,k),delta(ig,k),km(ig,k)
367         qpre=sqrt(q2(ig,k))
368         if (iflag_pbl.eq.8 ) then
369            if (aa(ig,k).gt.0.) then
370               q2(ig,k)=(qpre+aa(ig,k)*qpre*qpre)**2
371            else
372               q2(ig,k)=(qpre/(1.-aa(ig,k)*qpre))**2
373            endif
374         else ! iflag_pbl=9
375            if (aa(ig,k)*qpre.gt.0.9) then
376               q2(ig,k)=(qpre*10.)**2
377            else
378               q2(ig,k)=(qpre/(1.-aa(ig,k)*qpre))**2
379            endif
380         endif
381         q2(ig,k)=min(max(q2(ig,k),1.e-10),1.e4)
382c     print*,'Q2 L=',k,q2(ig,k),qpre*qpre
383                                                          enddo
384      enddo
385
386      endif ! Fin du cas 8
387
388c     print*,'OK8'
389
390c====================================================================
[766]391c   Calcul des coefficients de m�ange
[541]392c====================================================================
393      do k=2,klev
394c     print*,'k=',k
395                                                          do ig=1,ngrid
396cabde      print*,'KML=',l(ig,k),q2(ig,k),sm(ig,k)
397         zq=sqrt(q2(ig,k))
398         km(ig,k)=l(ig,k)*zq*sm(ig,k)
399         kn(ig,k)=km(ig,k)*alpha(ig,k)
400         kq(ig,k)=l(ig,k)*zq*0.2
401c     print*,'KML=',km(ig,k),kn(ig,k)
402                                                          enddo
403      enddo
404
[1403]405! Transport diffusif vertical de la TKE.
406      if (iflag_pbl.ge.9) then
407!       print*,'YAMADA VDIF'
408        q2(:,1)=q2(:,2)
409        call vdif_q2(dt,g,rconst,ngrid,plev,temp,kq,q2)
410      endif
[541]411
412c   Traitement des cas noctrunes avec l'introduction d'une longueur
413c   minilale.
414
415c====================================================================
416c   Traitement particulier pour les cas tres stables.
417c   D'apres Holtslag Boville.
418
[938]419      if (prt_level>1) THEN
420       print*,'YAMADA4 0'
421      endif !(prt_level>1) THEN
[541]422                                                          do ig=1,ngrid
423      coriol(ig)=1.e-4
424      pblhmin(ig)=0.07*ustar(ig)/max(abs(coriol(ig)),2.546e-5)
425                                                          enddo
426
[878]427!      print*,'pblhmin ',pblhmin
[541]428CTest a remettre 21 11 02
429c test abd 13 05 02      if(0.eq.1) then
430      if(1.eq.1) then
431      do k=2,klev
[789]432         do ig=1,ngrid
[541]433            if (teta(ig,2).gt.teta(ig,1)) then
434               qmin=ustar(ig)*(max(1.-zlev(ig,k)/pblhmin(ig),0.))**2
435               kmin=kap*zlev(ig,k)*qmin
436            else
437               kmin=-1. ! kmin n'est utilise que pour les SL stables.
438            endif
439            if (kn(ig,k).lt.kmin.or.km(ig,k).lt.kmin) then
440c               print*,'Seuil min Km K=',k,kmin,km(ig,k),kn(ig,k)
441c     s           ,sqrt(q2(ig,k)),pblhmin(ig),qmin/sm(ig,k)
442               kn(ig,k)=kmin
443               km(ig,k)=kmin
444               kq(ig,k)=kmin
445c   la longueur de melange est suposee etre l= kap z
446c   K=l q Sm d'ou q2=(K/l Sm)**2
447               q2(ig,k)=(qmin/sm(ig,k))**2
448            endif
449         enddo
450      enddo
451      endif
452
[938]453      if (prt_level>1) THEN
454       print*,'YAMADA4 1'
455      endif !(prt_level>1) THEN
[541]456c   Diagnostique pour stokage
457
[878]458      if(1.eq.0)then
[541]459      rino=rif
[878]460      smyam(1:ngrid,1)=0.
461      styam(1:ngrid,1)=0.
462      lyam(1:ngrid,1)=0.
463      knyam(1:ngrid,1)=0.
464      w2yam(1:ngrid,1)=0.
465      t2yam(1:ngrid,1)=0.
[541]466
[878]467      smyam(1:ngrid,2:klev)=sm(1:ngrid,2:klev)
468      styam(1:ngrid,2:klev)=sm(1:ngrid,2:klev)*alpha(1:ngrid,2:klev)
469      lyam(1:ngrid,2:klev)=l(1:ngrid,2:klev)
470      knyam(1:ngrid,2:klev)=kn(1:ngrid,2:klev)
471
[541]472c   Estimations de w'2 et T'2 d'apres Abdela et McFarlane
473
[878]474      w2yam(1:ngrid,2:klev)=q2(1:ngrid,2:klev)*0.24
475     s    +lyam(1:ngrid,2:klev)*5.17*kn(1:ngrid,2:klev)
476     s    *n2(1:ngrid,2:klev)/sqrt(q2(1:ngrid,2:klev))
[541]477
[878]478      t2yam(1:ngrid,2:klev)=9.1*kn(1:ngrid,2:klev)
479     s    *dtetadz(1:ngrid,2:klev)**2
480     s    /sqrt(q2(1:ngrid,2:klev))*lyam(1:ngrid,2:klev)
481      endif
[541]482
483c     print*,'OKFIN'
484      first=.false.
485      return
486      end
[1403]487      SUBROUTINE vdif_q2(timestep,gravity,rconst,ngrid,plev,temp,
488     &  kmy,q2)
489      use dimphy
490      IMPLICIT NONE
491c.......................................................................
492#include "dimensions.h"
493cccc#include "dimphy.h"
494c.......................................................................
495c
496c dt : pas de temps
497
498      real plev(klon,klev+1)
499      real temp(klon,klev)
500      real timestep
501      real gravity,rconst
502      real kstar(klon,klev+1),zz
503      real kmy(klon,klev+1)
504      real q2(klon,klev+1)
505      real deltap(klon,klev+1)
506      real denom(klon,klev+1),alpha(klon,klev+1),beta(klon,klev+1)
507      integer ngrid
508
509      integer i,k
510
511!       print*,'RD=',rconst
512      do k=1,klev
513         do i=1,ngrid
514c test
515!       print*,'i,k',i,k
516!       print*,'temp(i,k)=',temp(i,k)
517!       print*,'(plev(i,k)-plev(i,k+1))=',plev(i,k),plev(i,k+1)
518            zz=(plev(i,k)+plev(i,k+1))*gravity/(rconst*temp(i,k))
519            kstar(i,k)=0.125*(kmy(i,k+1)+kmy(i,k))*zz*zz
520     s      /(plev(i,k)-plev(i,k+1))*timestep
521         enddo
522      enddo
523
524      do k=2,klev
525         do i=1,ngrid
526            deltap(i,k)=0.5*(plev(i,k-1)-plev(i,k+1))
527         enddo
528      enddo
529      do i=1,ngrid
530         deltap(i,1)=0.5*(plev(i,1)-plev(i,2))
531         deltap(i,klev+1)=0.5*(plev(i,klev)-plev(i,klev+1))
532         denom(i,klev+1)=deltap(i,klev+1)+kstar(i,klev)
533         alpha(i,klev+1)=deltap(i,klev+1)*q2(i,klev+1)/denom(i,klev+1)
534         beta(i,klev+1)=kstar(i,klev)/denom(i,klev+1)
535      enddo
536
537      do k=klev,2,-1
538         do i=1,ngrid
539            denom(i,k)=deltap(i,k)+(1.-beta(i,k+1))*
540     s      kstar(i,k)+kstar(i,k-1)
541c   correction d'un bug 10 01 2001
542            alpha(i,k)=(q2(i,k)*deltap(i,k)
543     s      +kstar(i,k)*alpha(i,k+1))/denom(i,k)
544            beta(i,k)=kstar(i,k-1)/denom(i,k)
545         enddo
546      enddo
547
548c  Si on recalcule q2(1)
549      if(1.eq.0) then
550      do i=1,ngrid
551         denom(i,1)=deltap(i,1)+(1-beta(i,2))*kstar(i,1)
552         q2(i,1)=(q2(i,1)*deltap(i,1)
553     s      +kstar(i,1)*alpha(i,2))/denom(i,1)
554      enddo
555      endif
556c   sinon, on peut sauter cette boucle...
557
558      do k=2,klev+1
559         do i=1,ngrid
560            q2(i,k)=alpha(i,k)+beta(i,k)*q2(i,k-1)
561         enddo
562      enddo
563
564      return
565      end
566      SUBROUTINE vdif_q2e(timestep,gravity,rconst,ngrid,
567     &   plev,temp,kmy,q2)
568      use dimphy
569      IMPLICIT NONE
570c.......................................................................
571#include "dimensions.h"
572cccc#include "dimphy.h"
573c.......................................................................
574c
575c dt : pas de temps
576
577      real plev(klon,klev+1)
578      real temp(klon,klev)
579      real timestep
580      real gravity,rconst
581      real kstar(klon,klev+1),zz
582      real kmy(klon,klev+1)
583      real q2(klon,klev+1)
584      real deltap(klon,klev+1)
585      real denom(klon,klev+1),alpha(klon,klev+1),beta(klon,klev+1)
586      integer ngrid
587
588      integer i,k
589
590      do k=1,klev
591         do i=1,ngrid
592            zz=(plev(i,k)+plev(i,k+1))*gravity/(rconst*temp(i,k))
593            kstar(i,k)=0.125*(kmy(i,k+1)+kmy(i,k))*zz*zz
594     s      /(plev(i,k)-plev(i,k+1))*timestep
595         enddo
596      enddo
597
598      do k=2,klev
599         do i=1,ngrid
600            deltap(i,k)=0.5*(plev(i,k-1)-plev(i,k+1))
601         enddo
602      enddo
603      do i=1,ngrid
604         deltap(i,1)=0.5*(plev(i,1)-plev(i,2))
605         deltap(i,klev+1)=0.5*(plev(i,klev)-plev(i,klev+1))
606      enddo
607
608      do k=klev,2,-1
609         do i=1,ngrid
610            q2(i,k)=q2(i,k)+
611     s      ( kstar(i,k)*(q2(i,k+1)-q2(i,k))
612     s       -kstar(i,k-1)*(q2(i,k)-q2(i,k-1)) )
613     s      /deltap(i,k)
614         enddo
615      enddo
616
617      do i=1,ngrid
618         q2(i,1)=q2(i,1)+
619     s   ( kstar(i,1)*(q2(i,2)-q2(i,1))
620     s                                      )
621     s   /deltap(i,1)
622         q2(i,klev+1)=q2(i,klev+1)+
623     s   (
624     s    -kstar(i,klev)*(q2(i,klev+1)-q2(i,klev)) )
625     s   /deltap(i,klev+1)
626      enddo
627
628      return
629      end
Note: See TracBrowser for help on using the repository browser.