source: LMDZ4/trunk/libf/phylmd/orografi.F @ 524

Last change on this file since 524 was 524, checked in by lmdzadmin, 20 years ago

Initial revision

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 47.4 KB
Line 
1!
2! $Header$
3!
4      SUBROUTINE drag_noro (nlon,nlev,dtime,paprs,pplay,
5     e                   pmea,pstd, psig, pgam, pthe,ppic,pval,
6     e                   kgwd,kdx,ktest,
7     e                   t, u, v,
8     s                   pulow, pvlow, pustr, pvstr,
9     s                   d_t, d_u, d_v)
10c
11      IMPLICIT none
12c======================================================================
13c Auteur(s): F.Lott (LMD/CNRS) date: 19950201
14c Objet: Frottement de la montagne Interface
15c======================================================================
16c Arguments:
17c dtime---input-R- pas d'integration (s)
18c paprs---input-R-pression pour chaque inter-couche (en Pa)
19c pplay---input-R-pression pour le mileu de chaque couche (en Pa)
20c t-------input-R-temperature (K)
21c u-------input-R-vitesse horizontale (m/s)
22c v-------input-R-vitesse horizontale (m/s)
23c
24c d_t-----output-R-increment de la temperature             
25c d_u-----output-R-increment de la vitesse u
26c d_v-----output-R-increment de la vitesse v
27c======================================================================
28#include "dimensions.h"
29#include "dimphy.h"
30#include "YOMCST.h"
31c
32c ARGUMENTS
33c
34      INTEGER nlon,nlev
35      REAL dtime
36      REAL paprs(klon,klev+1)
37      REAL pplay(klon,klev)
38      REAL pmea(nlon),pstd(nlon),psig(nlon),pgam(nlon),pthe(nlon)
39      REAL ppic(nlon),pval(nlon)
40      REAL pulow(nlon),pvlow(nlon),pustr(nlon),pvstr(nlon)
41      REAL t(nlon,nlev), u(nlon,nlev), v(nlon,nlev)
42      REAL d_t(nlon,nlev), d_u(nlon,nlev), d_v(nlon,nlev)
43c
44      INTEGER i, k, kgwd, kdx(nlon), ktest(nlon)
45c
46c Variables locales:
47c
48      REAL zgeom(klon,klev)
49      REAL pdtdt(klon,klev), pdudt(klon,klev), pdvdt(klon,klev)
50      REAL pt(klon,klev), pu(klon,klev), pv(klon,klev)
51      REAL papmf(klon,klev),papmh(klon,klev+1)
52c
53c initialiser les variables de sortie (pour securite)
54c
55      DO i = 1,klon
56         pulow(i) = 0.0
57         pvlow(i) = 0.0
58         pustr(i) = 0.0
59         pvstr(i) = 0.0
60      ENDDO
61      DO k = 1, klev
62      DO i = 1, klon
63         d_t(i,k) = 0.0
64         d_u(i,k) = 0.0
65         d_v(i,k) = 0.0
66         pdudt(i,k)=0.0
67         pdvdt(i,k)=0.0
68         pdtdt(i,k)=0.0
69      ENDDO
70      ENDDO
71c
72c preparer les variables d'entree (attention: l'ordre des niveaux
73c verticaux augmente du haut vers le bas)
74c
75      DO k = 1, klev
76      DO i = 1, klon
77         pt(i,k) = t(i,klev-k+1)
78         pu(i,k) = u(i,klev-k+1)
79         pv(i,k) = v(i,klev-k+1)
80         papmf(i,k) = pplay(i,klev-k+1)
81      ENDDO
82      ENDDO
83      DO k = 1, klev+1
84      DO i = 1, klon
85         papmh(i,k) = paprs(i,klev-k+2)
86      ENDDO
87      ENDDO
88      DO i = 1, klon
89         zgeom(i,klev) = RD * pt(i,klev)
90     .                  * LOG(papmh(i,klev+1)/papmf(i,klev))
91      ENDDO
92      DO k = klev-1, 1, -1
93      DO i = 1, klon
94         zgeom(i,k) = zgeom(i,k+1) + RD * (pt(i,k)+pt(i,k+1))/2.0
95     .               * LOG(papmf(i,k+1)/papmf(i,k))
96      ENDDO
97      ENDDO
98c
99c appeler la routine principale
100c
101      CALL orodrag(klon,klev,kgwd,kdx,ktest,
102     .            dtime,
103     .            papmh, papmf, zgeom,
104     .            pt, pu, pv,
105     .            pmea, pstd, psig, pgam, pthe, ppic,pval,
106     .            pulow,pvlow,
107     .            pdudt,pdvdt,pdtdt)
108C
109      DO k = 1, klev
110      DO i = 1, klon
111         d_u(i,klev+1-k) = dtime*pdudt(i,k)
112         d_v(i,klev+1-k) = dtime*pdvdt(i,k)
113         d_t(i,klev+1-k) = dtime*pdtdt(i,k)
114         pustr(i)        = pustr(i)
115     .                    +rg*pdudt(i,k)*(papmh(i,k+1)-papmh(i,k))
116         pvstr(i)        = pvstr(i)
117     .                    +rg*pdvdt(i,k)*(papmh(i,k+1)-papmh(i,k))
118      ENDDO
119      ENDDO
120c
121      RETURN
122      END
123      SUBROUTINE orodrag( nlon,nlev
124     i                 , kgwd, kdx, ktest
125     r                 , ptsphy
126     r                 , paphm1,papm1,pgeom1,ptm1,pum1,pvm1
127     r                 , pmea, pstd, psig, pgamma, ptheta, ppic, pval
128c outputs
129     r                 , pulow,pvlow
130     r                 , pvom,pvol,pte )
131
132      implicit none
133
134c
135c
136c**** *gwdrag* - does the gravity wave parametrization.
137c
138c     purpose.
139c     --------
140c
141c          this routine computes the physical tendencies of the
142c     prognostic variables u,v  and t due to  vertical transports by
143c     subgridscale orographically excited gravity waves
144c
145c**   interface.
146c     ----------
147c          called from *callpar*.
148c
149c          the routine takes its input from the long-term storage:
150c          u,v,t and p at t-1.
151c
152c        explicit arguments :
153c        --------------------
154c     ==== inputs ===
155c     ==== outputs ===
156c
157c        implicit arguments :   none
158c        --------------------
159c
160c      implicit logical (l)
161c
162c     method.
163c     -------
164c
165c     externals.
166c     ----------
167      integer ismin, ismax
168      external ismin, ismax
169c
170c     reference.
171c     ----------
172c
173c     author.
174c     -------
175c     m.miller + b.ritter   e.c.m.w.f.     15/06/86.
176c
177c     f.lott + m. miller    e.c.m.w.f.     22/11/94
178c-----------------------------------------------------------------------
179c
180c
181#include "dimensions.h"
182#include "dimphy.h"
183#include "YOMCST.h"
184#include "YOEGWD.h"
185c-----------------------------------------------------------------------
186c
187c*       0.1   arguments
188c              ---------
189c
190c
191      integer nlon, nlev, klevm1
192      integer kgwd, jl, ilevp1, jk, ji
193      real zdelp, ztemp, zforc, ztend
194      real rover, zb, zc, zconb, zabsv
195      real zzd1, ratio, zbet, zust,zvst, zdis
196      real  pte(nlon,nlev),
197     *      pvol(nlon,nlev),
198     *      pvom(nlon,nlev),
199     *      pulow(klon),
200     *      pvlow(klon)
201      real  pum1(nlon,nlev),
202     *      pvm1(nlon,nlev),
203     *      ptm1(nlon,nlev),
204     *      pmea(nlon),pstd(nlon),psig(nlon),
205     *      pgamma(nlon),ptheta(nlon),ppic(nlon),pval(nlon),
206     *      pgeom1(nlon,nlev),
207     *      papm1(nlon,nlev),
208     *      paphm1(nlon,nlev+1)
209c
210      integer  kdx(nlon),ktest(nlon)
211c-----------------------------------------------------------------------
212c
213c*       0.2   local arrays
214c              ------------
215      integer  isect(klon),
216     *         icrit(klon),
217     *         ikcrith(klon),
218     *         ikenvh(klon),
219     *         iknu(klon),
220     *         iknu2(klon),
221     *         ikcrit(klon),
222     *         ikhlim(klon)
223c
224      real   ztau(klon,klev+1),
225     $       ztauf(klon,klev+1),
226     *       zstab(klon,klev+1),
227     *       zvph(klon,klev+1),
228     *       zrho(klon,klev+1),
229     *       zri(klon,klev+1),
230     *       zpsi(klon,klev+1),
231     *       zzdep(klon,klev)
232      real   zdudt(klon),
233     *       zdvdt(klon),
234     *       zdtdt(klon),
235     *       zdedt(klon),
236     *       zvidis(klon),
237     *       znu(klon),
238     *       zd1(klon),
239     *       zd2(klon),
240     *       zdmod(klon)
241      real ztmst, ptsphy, zrtmst
242c
243c------------------------------------------------------------------
244c
245c*         1.    initialization
246c                --------------
247c
248 100  continue
249c
250c     ------------------------------------------------------------------
251c
252c*         1.1   computational constants
253c                -----------------------
254c
255 110  continue
256c
257c     ztmst=twodt
258c     if(nstep.eq.nstart) ztmst=0.5*twodt
259      klevm1=klev-1
260      ztmst=ptsphy
261      zrtmst=1./ztmst
262c     ------------------------------------------------------------------
263c
264 120  continue
265c
266c     ------------------------------------------------------------------
267c
268c*         1.3   check whether row contains point for printing
269c                ---------------------------------------------
270c
271 130  continue
272c
273c     ------------------------------------------------------------------
274c
275c*         2.     precompute basic state variables.
276c*                ---------- ----- ----- ----------
277c*                define low level wind, project winds in plane of
278c*                low level wind, determine sector in which to take
279c*                the variance and set indicator for critical levels.
280c
281  200 continue
282c
283c
284c
285      call orosetup
286     *     ( nlon, ktest
287     *     , ikcrit, ikcrith, icrit,  ikenvh,iknu,iknu2
288     *     , paphm1, papm1 , pum1   , pvm1 , ptm1 , pgeom1, pstd
289     *     , zrho  , zri   , zstab  , ztau , zvph , zpsi, zzdep
290     *     , pulow, pvlow
291     *     , ptheta,pgamma,pmea,ppic,pval,znu  ,zd1,  zd2,  zdmod )
292c
293c
294c
295c***********************************************************
296c
297c
298c*         3.      compute low level stresses using subcritical and
299c*                 supercritical forms.computes anisotropy coefficient
300c*                 as measure of orographic twodimensionality.
301c
302  300 continue
303c
304      call gwstress
305     *    ( nlon  , nlev
306     *    , ktest , icrit, ikenvh, iknu
307     *    , zrho  , zstab, zvph  , pstd,  psig, pmea, ppic
308     *    , ztau
309     *    , pgeom1,zdmod)
310c
311c
312c*         4.      compute stress profile.
313c*                 ------- ------ --------
314c
315  400 continue
316c
317c
318      call gwprofil
319     *       (  nlon , nlev
320     *       , kgwd   , kdx , ktest
321     *       , ikcrith, icrit
322     *       , paphm1, zrho   , zstab ,  zvph
323     *       , zri   , ztau   
324     *       , zdmod , psig  , pstd)
325c
326c
327c*         5.      compute tendencies.
328c*                 -------------------
329c
330  500 continue
331c
332c  explicit solution at all levels for the gravity wave
333c  implicit solution for the blocked levels
334
335      do 510 jl=kidia,kfdia
336      zvidis(jl)=0.0
337      zdudt(jl)=0.0
338      zdvdt(jl)=0.0
339      zdtdt(jl)=0.0
340  510 continue
341c
342      ilevp1=klev+1
343c
344c
345      do 524 jk=1,klev
346c
347c
348c     do 523 jl=1,kgwd
349c     ji=kdx(jl)
350c  Modif vectorisation 02/04/2004
351      do 523 ji=kidia,kfdia
352      if(ktest(ji).eq.1) then
353
354      zdelp=paphm1(ji,jk+1)-paphm1(ji,jk)
355      ztemp=-rg*(ztau(ji,jk+1)-ztau(ji,jk))/(zvph(ji,ilevp1)*zdelp)
356      zdudt(ji)=(pulow(ji)*zd1(ji)-pvlow(ji)*zd2(ji))*ztemp/zdmod(ji)
357      zdvdt(ji)=(pvlow(ji)*zd1(ji)+pulow(ji)*zd2(ji))*ztemp/zdmod(ji)
358c
359c controle des overshoots:
360c
361      zforc=sqrt(zdudt(ji)**2+zdvdt(ji)**2)+1.E-12
362      ztend=sqrt(pum1(ji,jk)**2+pvm1(ji,jk)**2)/ztmst+1.E-12
363      rover=0.25
364      if(zforc.ge.rover*ztend)then
365        zdudt(ji)=rover*ztend/zforc*zdudt(ji)
366        zdvdt(ji)=rover*ztend/zforc*zdvdt(ji)
367      endif
368c
369c fin du controle des overshoots
370c
371      if(jk.ge.ikenvh(ji)) then
372         zb=1.0-0.18*pgamma(ji)-0.04*pgamma(ji)**2
373         zc=0.48*pgamma(ji)+0.3*pgamma(ji)**2
374         zconb=2.*ztmst*gkwake*psig(ji)/(4.*pstd(ji))
375         zabsv=sqrt(pum1(ji,jk)**2+pvm1(ji,jk)**2)/2.
376         zzd1=zb*cos(zpsi(ji,jk))**2+zc*sin(zpsi(ji,jk))**2
377             ratio=(cos(zpsi(ji,jk))**2+pgamma(ji)*sin(zpsi(ji,jk))**2)/
378     *   (pgamma(ji)*cos(zpsi(ji,jk))**2+sin(zpsi(ji,jk))**2)
379         zbet=max(0.,2.-1./ratio)*zconb*zzdep(ji,jk)*zzd1*zabsv
380c
381c simplement oppose au vent
382c
383         zdudt(ji)=-pum1(ji,jk)/ztmst
384         zdvdt(ji)=-pvm1(ji,jk)/ztmst
385c
386c  projection dans la direction de l'axe principal de l'orographie
387cmod     zdudt(ji)=-(pum1(ji,jk)*cos(ptheta(ji)*rpi/180.)
388cmod *              +pvm1(ji,jk)*sin(ptheta(ji)*rpi/180.))
389cmod *              *cos(ptheta(ji)*rpi/180.)/ztmst
390cmod     zdvdt(ji)=-(pum1(ji,jk)*cos(ptheta(ji)*rpi/180.)
391cmod *              +pvm1(ji,jk)*sin(ptheta(ji)*rpi/180.))
392cmod *              *sin(ptheta(ji)*rpi/180.)/ztmst
393         zdudt(ji)=zdudt(ji)*(zbet/(1.+zbet))
394         zdvdt(ji)=zdvdt(ji)*(zbet/(1.+zbet))
395      end if
396      pvom(ji,jk)=zdudt(ji)
397      pvol(ji,jk)=zdvdt(ji)
398      zust=pum1(ji,jk)+ztmst*zdudt(ji)
399      zvst=pvm1(ji,jk)+ztmst*zdvdt(ji)
400      zdis=0.5*(pum1(ji,jk)**2+pvm1(ji,jk)**2-zust**2-zvst**2)
401      zdedt(ji)=zdis/ztmst
402      zvidis(ji)=zvidis(ji)+zdis*zdelp
403      zdtdt(ji)=zdedt(ji)/rcpd
404c     pte(ji,jk)=zdtdt(ji)
405c
406c  ENCORE UN TRUC POUR EVITER LES EXPLOSIONS
407c
408      pte(ji,jk)=0.0
409
410      endif
411  523 continue
412
413  524 continue
414c
415c
416      return
417      end
418      SUBROUTINE orosetup
419     *         ( nlon   , ktest
420     *         , kkcrit, kkcrith, kcrit
421     *         , kkenvh, kknu  , kknu2
422     *         , paphm1, papm1 , pum1   , pvm1 , ptm1  , pgeom1, pstd
423     *         , prho  , pri   , pstab  , ptau , pvph  ,ppsi, pzdep
424     *         , pulow , pvlow 
425     *         , ptheta, pgamma, pmea, ppic, pval
426     *         , pnu  ,  pd1  ,  pd2  ,pdmod  )
427c
428c**** *gwsetup*
429c
430c     purpose.
431c     --------
432c
433c**   interface.
434c     ----------
435c          from *orodrag*
436c
437c        explicit arguments :
438c        --------------------
439c     ==== inputs ===
440c     ==== outputs ===
441c
442c        implicit arguments :   none
443c        --------------------
444c
445c     method.
446c     -------
447c
448c
449c     externals.
450c     ----------
451c
452c
453c     reference.
454c     ----------
455c
456c        see ecmwf research department documentation of the "i.f.s."
457c
458c     author.
459c     -------
460c
461c     modifications.
462c     --------------
463c     f.lott  for the new-gwdrag scheme november 1993
464c
465c-----------------------------------------------------------------------
466      implicit none
467c
468
469#include "dimensions.h"
470#include "dimphy.h"
471#include "YOMCST.h"
472#include "YOEGWD.h"
473
474c-----------------------------------------------------------------------
475c
476c*       0.1   arguments
477c              ---------
478c
479      integer nlon
480      integer jl, jk
481      real zdelp
482
483      integer kkcrit(nlon),kkcrith(nlon),kcrit(nlon),
484     *        ktest(nlon),kkenvh(nlon)
485
486c
487      real paphm1(nlon,klev+1),papm1(nlon,klev),pum1(nlon,klev),
488     *     pvm1(nlon,klev),ptm1(nlon,klev),pgeom1(nlon,klev),
489     *     prho(nlon,klev+1),pri(nlon,klev+1),pstab(nlon,klev+1),
490     *     ptau(nlon,klev+1),pvph(nlon,klev+1),ppsi(nlon,klev+1),
491     *     pzdep(nlon,klev)
492       real pulow(nlon),pvlow(nlon),ptheta(nlon),pgamma(nlon),pnu(nlon),
493     *     pd1(nlon),pd2(nlon),pdmod(nlon)
494      real pstd(nlon),pmea(nlon),ppic(nlon),pval(nlon)
495c
496c-----------------------------------------------------------------------
497c
498c*       0.2   local arrays
499c              ------------
500c
501c
502      integer ilevm1, ilevm2, ilevh
503      real zcons1, zcons2,zcons3, zhgeo
504      real zu, zphi, zvt1,zvt2, zst, zvar, zdwind, zwind
505      real zstabm, zstabp, zrhom,  zrhop, alpha
506      real zggeenv, zggeom1,zgvar
507      logical lo
508      logical ll1(klon,klev+1)
509      integer kknu(klon),kknu2(klon),kknub(klon),kknul(klon),
510     *        kentp(klon),ncount(klon) 
511c
512      real zhcrit(klon,klev),zvpf(klon,klev),
513     *     zdp(klon,klev)
514      real znorm(klon),zb(klon),zc(klon),
515     *      zulow(klon),zvlow(klon),znup(klon),znum(klon)
516c
517c     ------------------------------------------------------------------
518c
519c*         1.    initialization
520c                --------------
521c
522c     print *,' entree gwsetup'
523 100  continue
524c
525c     ------------------------------------------------------------------
526c
527c*         1.1   computational constants
528c                -----------------------
529c
530 110  continue
531c
532      ilevm1=klev-1
533      ilevm2=klev-2
534      ilevh =klev/3
535c
536      zcons1=1./rd
537cold  zcons2=g**2/cpd
538      zcons2=rg**2/rcpd
539cold  zcons3=1.5*api
540      zcons3=1.5*rpi
541c
542c
543c     ------------------------------------------------------------------
544c
545c*         2.
546c                --------------
547c
548 200  continue
549c
550c     ------------------------------------------------------------------
551c
552c*         2.1     define low level wind, project winds in plane of
553c*                 low level wind, determine sector in which to take
554c*                 the variance and set indicator for critical levels.
555c
556c
557c
558      do 2001 jl=kidia,kfdia
559      kknu(jl)    =klev
560      kknu2(jl)   =klev
561      kknub(jl)   =klev
562      kknul(jl)   =klev
563      pgamma(jl) =max(pgamma(jl),gtsec)
564      ll1(jl,klev+1)=.false.
565 2001 continue
566c
567c Ajouter une initialisation (L. Li, le 23fev99):
568c
569      do jk=klev,ilevh,-1
570      do jl=kidia,kfdia
571      ll1(jl,jk)= .FALSE.
572      ENDDO
573      ENDDO
574c
575c*      define top of low level flow
576c       ----------------------------
577      do 2002 jk=klev,ilevh,-1
578      do 2003 jl=kidia,kfdia
579      lo=(paphm1(jl,jk)/paphm1(jl,klev+1)).ge.gsigcr
580      if(lo) then
581        kkcrit(jl)=jk
582      endif
583      zhcrit(jl,jk)=ppic(jl)
584      zhgeo=pgeom1(jl,jk)/rg
585      ll1(jl,jk)=(zhgeo.gt.zhcrit(jl,jk))
586      if(ll1(jl,jk).xor.ll1(jl,jk+1)) then
587        kknu(jl)=jk
588      endif
589      if(.not.ll1(jl,ilevh))kknu(jl)=ilevh
590 2003 continue
591 2002 continue
592      do 2004 jk=klev,ilevh,-1
593      do 2005 jl=kidia,kfdia
594      zhcrit(jl,jk)=ppic(jl)-pval(jl)
595      zhgeo=pgeom1(jl,jk)/rg
596      ll1(jl,jk)=(zhgeo.gt.zhcrit(jl,jk))
597      if(ll1(jl,jk).xor.ll1(jl,jk+1)) then
598        kknu2(jl)=jk
599      endif
600      if(.not.ll1(jl,ilevh))kknu2(jl)=ilevh
601 2005 continue
602 2004 continue
603      do 2006 jk=klev,ilevh,-1
604      do 2007 jl=kidia,kfdia
605      zhcrit(jl,jk)=amax1(ppic(jl)-pmea(jl),pmea(jl)-pval(jl))
606      zhgeo=pgeom1(jl,jk)/rg
607      ll1(jl,jk)=(zhgeo.gt.zhcrit(jl,jk))
608      if(ll1(jl,jk).xor.ll1(jl,jk+1)) then
609        kknub(jl)=jk
610      endif
611      if(.not.ll1(jl,ilevh))kknub(jl)=ilevh
612 2007 continue
613 2006 continue
614c
615      do 2010 jl=kidia,kfdia 
616      kknu(jl)=min(kknu(jl),nktopg)
617      kknu2(jl)=min(kknu2(jl),nktopg)
618      kknub(jl)=min(kknub(jl),nktopg)
619      kknul(jl)=klev
620 2010 continue     
621c
622
623 210  continue
624c
625c
626cc*     initialize various arrays
627c
628      do 2107 jl=kidia,kfdia
629      prho(jl,klev+1)  =0.0
630      pstab(jl,klev+1) =0.0
631      pstab(jl,1)      =0.0
632      pri(jl,klev+1)   =9999.0
633      ppsi(jl,klev+1)  =0.0
634      pri(jl,1)        =0.0
635      pvph(jl,1)       =0.0
636      pulow(jl)        =0.0
637      pvlow(jl)        =0.0
638      zulow(jl)        =0.0
639      zvlow(jl)        =0.0
640      kkcrith(jl)      =klev
641      kkenvh(jl)       =klev
642      kentp(jl)        =klev
643      kcrit(jl)        =1
644      ncount(jl)       =0
645      ll1(jl,klev+1)   =.false.
646 2107 continue
647c
648c*     define low-level flow
649c      ---------------------
650c
651      do 223 jk=klev,2,-1
652      do 222 jl=kidia,kfdia
653      if(ktest(jl).eq.1) then
654        zdp(jl,jk)=papm1(jl,jk)-papm1(jl,jk-1)
655        prho(jl,jk)=2.*paphm1(jl,jk)*zcons1/(ptm1(jl,jk)+ptm1(jl,jk-1))
656        pstab(jl,jk)=2.*zcons2/(ptm1(jl,jk)+ptm1(jl,jk-1))*
657     *  (1.-rcpd*prho(jl,jk)*(ptm1(jl,jk)-ptm1(jl,jk-1))/zdp(jl,jk))
658        pstab(jl,jk)=max(pstab(jl,jk),gssec)
659      endif
660  222 continue
661  223 continue
662c
663c********************************************************************
664c
665c*     define blocked flow
666c      -------------------
667      do 2115 jk=klev,ilevh,-1
668      do 2116 jl=kidia,kfdia
669      if(jk.ge.kknub(jl).and.jk.le.kknul(jl)) then
670        pulow(jl)=pulow(jl)+pum1(jl,jk)*(paphm1(jl,jk+1)-paphm1(jl,jk))
671        pvlow(jl)=pvlow(jl)+pvm1(jl,jk)*(paphm1(jl,jk+1)-paphm1(jl,jk))
672      end if
673 2116 continue
674 2115 continue
675      do 2110 jl=kidia,kfdia
676      pulow(jl)=pulow(jl)/(paphm1(jl,kknul(jl)+1)-paphm1(jl,kknub(jl)))
677      pvlow(jl)=pvlow(jl)/(paphm1(jl,kknul(jl)+1)-paphm1(jl,kknub(jl)))
678      znorm(jl)=max(sqrt(pulow(jl)**2+pvlow(jl)**2),gvsec)
679      pvph(jl,klev+1)=znorm(jl)
680 2110 continue
681c
682c*******  setup orography axes and define plane of profiles  *******
683c
684      do 2112 jl=kidia,kfdia
685      lo=(pulow(jl).lt.gvsec).and.(pulow(jl).ge.-gvsec)
686      if(lo) then
687        zu=pulow(jl)+2.*gvsec
688      else
689        zu=pulow(jl)
690      endif
691      zphi=atan(pvlow(jl)/zu)
692      ppsi(jl,klev+1)=ptheta(jl)*rpi/180.-zphi
693      zb(jl)=1.-0.18*pgamma(jl)-0.04*pgamma(jl)**2
694      zc(jl)=0.48*pgamma(jl)+0.3*pgamma(jl)**2
695      pd1(jl)=zb(jl)-(zb(jl)-zc(jl))*(sin(ppsi(jl,klev+1))**2)
696      pd2(jl)=(zb(jl)-zc(jl))*sin(ppsi(jl,klev+1))*cos(ppsi(jl,klev+1))
697      pdmod(jl)=sqrt(pd1(jl)**2+pd2(jl)**2)
698 2112 continue
699c
700c  ************ define flow in plane of lowlevel stress *************
701c
702      do 213 jk=1,klev
703      do 212 jl=kidia,kfdia
704      if(ktest(jl).eq.1)  then
705        zvt1       =pulow(jl)*pum1(jl,jk)+pvlow(jl)*pvm1(jl,jk)
706        zvt2       =-pvlow(jl)*pum1(jl,jk)+pulow(jl)*pvm1(jl,jk)
707        zvpf(jl,jk)=(zvt1*pd1(jl)+zvt2*pd2(jl))/(znorm(jl)*pdmod(jl))
708      endif
709      ptau(jl,jk)  =0.0
710      pzdep(jl,jk) =0.0
711      ppsi(jl,jk)  =0.0
712      ll1(jl,jk)   =.false.
713  212 continue
714  213 continue
715      do 215 jk=2,klev
716      do 214 jl=kidia,kfdia
717      if(ktest(jl).eq.1) then
718        zdp(jl,jk)=papm1(jl,jk)-papm1(jl,jk-1)
719        pvph(jl,jk)=((paphm1(jl,jk)-papm1(jl,jk-1))*zvpf(jl,jk)+
720     *            (papm1(jl,jk)-paphm1(jl,jk))*zvpf(jl,jk-1))
721     *            /zdp(jl,jk)
722        if(pvph(jl,jk).lt.gvsec) then
723          pvph(jl,jk)=gvsec
724          kcrit(jl)=jk
725        endif
726      endif
727  214 continue
728  215 continue
729c
730c
731c*         2.2     brunt-vaisala frequency and density at half levels.
732c
733  220 continue
734c
735      do 2211 jk=ilevh,klev
736      do 221 jl=kidia,kfdia
737      if(ktest(jl).eq.1) then
738      if(jk.ge.(kknub(jl)+1).and.jk.le.kknul(jl)) then
739           zst=zcons2/ptm1(jl,jk)*(1.-rcpd*prho(jl,jk)*
740     *                   (ptm1(jl,jk)-ptm1(jl,jk-1))/zdp(jl,jk))
741           pstab(jl,klev+1)=pstab(jl,klev+1)+zst*zdp(jl,jk)
742           pstab(jl,klev+1)=max(pstab(jl,klev+1),gssec)
743           prho(jl,klev+1)=prho(jl,klev+1)+paphm1(jl,jk)*2.*zdp(jl,jk)
744     *                   *zcons1/(ptm1(jl,jk)+ptm1(jl,jk-1))
745      endif
746      endif
747  221 continue
748 2211 continue
749c
750      do 2212 jl=kidia,kfdia
751        pstab(jl,klev+1)=pstab(jl,klev+1)/(papm1(jl,kknul(jl))
752     *                                          -papm1(jl,kknub(jl)))
753        prho(jl,klev+1)=prho(jl,klev+1)/(papm1(jl,kknul(jl))
754     *                                          -papm1(jl,kknub(jl)))
755        zvar=pstd(jl)
756 2212 continue
757c
758c*         2.3     mean flow richardson number.
759c*                 and critical height for froude layer
760c
761  230 continue
762c
763      do 232 jk=2,klev
764      do 231 jl=kidia,kfdia
765      if(ktest(jl).eq.1) then
766        zdwind=max(abs(zvpf(jl,jk)-zvpf(jl,jk-1)),gvsec)
767        pri(jl,jk)=pstab(jl,jk)*(zdp(jl,jk)
768     *          /(rg*prho(jl,jk)*zdwind))**2
769        pri(jl,jk)=max(pri(jl,jk),grcrit)
770      endif
771  231 continue
772  232 continue
773 
774c
775c
776c*      define top of 'envelope' layer
777c       ----------------------------
778
779      do 233 jl=kidia,kfdia
780      pnu (jl)=0.0
781      znum(jl)=0.0
782 233  continue
783     
784      do 234 jk=2,klev-1
785      do 234 jl=kidia,kfdia
786     
787      if(ktest(jl).eq.1) then
788       
789      if (jk.ge.kknub(jl)) then
790         
791            znum(jl)=pnu(jl)
792            zwind=(pulow(jl)*pum1(jl,jk)+pvlow(jl)*pvm1(jl,jk))/
793     *            max(sqrt(pulow(jl)**2+pvlow(jl)**2),gvsec)
794            zwind=max(sqrt(zwind**2),gvsec)
795            zdelp=paphm1(jl,jk+1)-paphm1(jl,jk)
796            zstabm=sqrt(max(pstab(jl,jk  ),gssec))
797            zstabp=sqrt(max(pstab(jl,jk+1),gssec))
798            zrhom=prho(jl,jk  )
799            zrhop=prho(jl,jk+1)
800            pnu(jl) = pnu(jl) + (zdelp/rg)*
801     *            ((zstabp/zrhop+zstabm/zrhom)/2.)/zwind     
802            if((znum(jl).le.gfrcrit).and.(pnu(jl).gt.gfrcrit)
803     *                          .and.(kkenvh(jl).eq.klev))
804     *      kkenvh(jl)=jk
805     
806      endif   
807
808      endif
809     
810 234  continue
811     
812c  calculation of a dynamical mixing height for the breaking
813c  of gravity waves:
814
815             
816      do 235 jl=kidia,kfdia
817      znup(jl)=0.0
818      znum(jl)=0.0
819 235  continue
820
821      do 236 jk=klev-1,2,-1
822      do 236 jl=kidia,kfdia
823     
824      if(ktest(jl).eq.1) then
825
826            znum(jl)=znup(jl)
827            zwind=(pulow(jl)*pum1(jl,jk)+pvlow(jl)*pvm1(jl,jk))/
828     *            max(sqrt(pulow(jl)**2+pvlow(jl)**2),gvsec)
829            zwind=max(sqrt(zwind**2),gvsec)
830            zdelp=paphm1(jl,jk+1)-paphm1(jl,jk)
831            zstabm=sqrt(max(pstab(jl,jk  ),gssec))
832            zstabp=sqrt(max(pstab(jl,jk+1),gssec))
833            zrhom=prho(jl,jk  )
834            zrhop=prho(jl,jk+1)
835            znup(jl) = znup(jl) + (zdelp/rg)*
836     *            ((zstabp/zrhop+zstabm/zrhom)/2.)/zwind     
837            if((znum(jl).le.rpi/2.).and.(znup(jl).gt.rpi/2.)
838     *                          .and.(kkcrith(jl).eq.klev))
839     *      kkcrith(jl)=jk
840     
841      endif
842     
843 236  continue
844 
845      do 237 jl=kidia,kfdia
846      kkcrith(jl)=min0(kkcrith(jl),kknu2(jl))
847      kkcrith(jl)=max0(kkcrith(jl),ilevh*2)
848 237  continue         
849c
850c     directional info for flow blocking *************************
851c
852      do 251 jk=ilevh,klev   
853      do 252 jl=kidia,kfdia
854      if(jk.ge.kkenvh(jl)) then
855      lo=(pum1(jl,jk).lt.gvsec).and.(pum1(jl,jk).ge.-gvsec)
856      if(lo) then
857        zu=pum1(jl,jk)+2.*gvsec
858      else
859        zu=pum1(jl,jk)
860      endif
861       zphi=atan(pvm1(jl,jk)/zu)
862       ppsi(jl,jk)=ptheta(jl)*rpi/180.-zphi
863      end if
864 252  continue
865 251  continue
866c      forms the vertical 'leakiness' **************************
867
868      alpha=3.
869     
870      do 254  jk=ilevh,klev
871      do 253  jl=kidia,kfdia
872      if(jk.ge.kkenvh(jl)) then
873        zggeenv=amax1(1.,
874     *          (pgeom1(jl,kkenvh(jl))+pgeom1(jl,kkenvh(jl)-1))/2.)     
875        zggeom1=amax1(pgeom1(jl,jk),1.)
876        zgvar=amax1(pstd(jl)*rg,1.)     
877cmod    pzdep(jl,jk)=sqrt((zggeenv-zggeom1)/(zggeom1+zgvar))     
878        pzdep(jl,jk)=(pgeom1(jl,kkenvh(jl)-1)-pgeom1(jl,  jk))/
879     *               (pgeom1(jl,kkenvh(jl)-1)-pgeom1(jl,klev))
880      end if
881 253  continue
882 254  continue
883
884 260  continue
885
886      return
887      end
888      SUBROUTINE gwstress
889     *         (  nlon  , nlev
890     *         , ktest, kcrit, kkenvh
891     *         , kknu
892     *         , prho  , pstab , pvph  , pstd, psig
893     *         , pmea , ppic  , ptau 
894     *         , pgeom1 , pdmod )
895c
896c**** *gwstress*
897c
898c     purpose.
899c     --------
900c
901c**   interface.
902c     ----------
903c     call *gwstress*  from *gwdrag*
904c
905c        explicit arguments :
906c        --------------------
907c     ==== inputs ===
908c     ==== outputs ===
909c
910c        implicit arguments :   none
911c        --------------------
912c
913c     method.
914c     -------
915c
916c
917c     externals.
918c     ----------
919c
920c
921c     reference.
922c     ----------
923c
924c        see ecmwf research department documentation of the "i.f.s."
925c
926c     author.
927c     -------
928c
929c     modifications.
930c     --------------
931c     f. lott put the new gwd on ifs      22/11/93
932c
933c-----------------------------------------------------------------------
934      implicit none
935#include "dimensions.h"
936#include "dimphy.h"
937#include "YOMCST.h"
938#include "YOEGWD.h"
939
940c-----------------------------------------------------------------------
941c
942c*       0.1   arguments
943c              ---------
944c
945      integer nlon, nlev
946      integer kcrit(nlon),
947     *        ktest(nlon),kkenvh(nlon),kknu(nlon)
948c
949      real prho(nlon,nlev+1),pstab(nlon,nlev+1),ptau(nlon,nlev+1),
950     *     pvph(nlon,nlev+1),
951     *     pgeom1(nlon,nlev),pstd(nlon)
952c
953      real psig(nlon)
954      real pmea(nlon),ppic(nlon)
955      real pdmod(nlon)
956c
957c-----------------------------------------------------------------------
958c
959c*       0.2   local arrays
960c              ------------
961      integer jl
962      real zblock, zvar, zeff
963      logical lo
964c
965c-----------------------------------------------------------------------
966c
967c*       0.3   functions
968c              ---------
969c     ------------------------------------------------------------------
970c
971c*         1.    initialization
972c                --------------
973c
974 100  continue
975c
976c*         3.1     gravity wave stress.
977c
978  300 continue
979c
980c
981      do 301 jl=kidia,kfdia
982      if(ktest(jl).eq.1) then
983     
984c  effective mountain height above the blocked flow
985 
986         if(kkenvh(jl).eq.klev)then
987         zblock=0.0
988         else
989         zblock=(pgeom1(jl,kkenvh(jl))+pgeom1(jl,kkenvh(jl)+1))/2./rg         
990         endif
991     
992        zvar=ppic(jl)-pmea(jl)
993        zeff=amax1(0.,zvar-zblock)
994
995        ptau(jl,klev+1)=prho(jl,klev+1)*gkdrag*psig(jl)*zeff**2
996     *    /4./pstd(jl)*pvph(jl,klev+1)*pdmod(jl)*sqrt(pstab(jl,klev+1))
997
998c  too small value of stress or  low level flow include critical level
999c  or low level flow:  gravity wave stress nul.
1000               
1001        lo=(ptau(jl,klev+1).lt.gtsec).or.(kcrit(jl).ge.kknu(jl))
1002     *      .or.(pvph(jl,klev+1).lt.gvcrit)
1003c       if(lo) ptau(jl,klev+1)=0.0
1004     
1005      else
1006     
1007          ptau(jl,klev+1)=0.0
1008         
1009      endif
1010     
1011  301 continue
1012c
1013      return
1014      end
1015      SUBROUTINE GWPROFIL
1016     *         ( NLON, NLEV
1017     *         , kgwd, kdx , ktest
1018     *         , KKCRITH, KCRIT
1019     *         , PAPHM1, PRHO   , PSTAB  , PVPH , PRI , PTAU
1020     *         , pdmod   , psig , pvar)
1021
1022C**** *GWPROFIL*
1023C
1024C     PURPOSE.
1025C     --------
1026C
1027C**   INTERFACE.
1028C     ----------
1029C          FROM *GWDRAG*
1030C
1031C        EXPLICIT ARGUMENTS :
1032C        --------------------
1033C     ==== INPUTS ===
1034C     ==== OUTPUTS ===
1035C
1036C        IMPLICIT ARGUMENTS :   NONE
1037C        --------------------
1038C
1039C     METHOD:
1040C     -------
1041C     THE STRESS PROFILE FOR GRAVITY WAVES IS COMPUTED AS FOLLOWS:
1042C     IT IS CONSTANT (NO GWD) AT THE LEVELS BETWEEN THE GROUND
1043C     AND THE TOP OF THE BLOCKED LAYER (KKENVH).
1044C     IT DECREASES LINEARLY WITH HEIGHTS FROM THE TOP OF THE
1045C     BLOCKED LAYER TO 3*VAROR (kKNU), TO SIMULATES LEE WAVES OR
1046C     NONLINEAR GRAVITY WAVE BREAKING.
1047C     ABOVE IT IS CONSTANT, EXCEPT WHEN THE WAVE ENCOUNTERS A CRITICAL
1048C     LEVEL (KCRIT) OR WHEN IT BREAKS.
1049C     
1050C
1051C
1052C     EXTERNALS.
1053C     ----------
1054C
1055C
1056C     REFERENCE.
1057C     ----------
1058C
1059C        SEE ECMWF RESEARCH DEPARTMENT DOCUMENTATION OF THE "I.F.S."
1060C
1061C     AUTHOR.
1062C     -------
1063C
1064C     MODIFICATIONS.
1065C     --------------
1066C     PASSAGE OF THE NEW GWDRAG TO I.F.S. (F. LOTT, 22/11/93)
1067C-----------------------------------------------------------------------
1068      implicit none
1069C
1070
1071C
1072
1073#include "dimensions.h"
1074#include "dimphy.h"
1075#include "YOMCST.h"
1076#include "YOEGWD.h"
1077
1078C-----------------------------------------------------------------------
1079C
1080C*       0.1   ARGUMENTS
1081C              ---------
1082C
1083      integer nlon,nlev
1084      INTEGER KKCRITH(NLON),KCRIT(NLON)
1085     *       ,kdx(nlon) , ktest(nlon)
1086
1087C
1088      REAL PAPHM1(NLON,NLEV+1), PSTAB(NLON,NLEV+1),
1089     *     PRHO  (NLON,NLEV+1), PVPH (NLON,NLEV+1),
1090     *     PRI   (NLON,NLEV+1), PTAU(NLON,NLEV+1)
1091     
1092      REAL pdmod (NLON) , psig(NLON),
1093     *     pvar(NLON)
1094     
1095C-----------------------------------------------------------------------
1096C
1097C*       0.2   LOCAL ARRAYS
1098C              ------------
1099C
1100      integer ilevh, ji, kgwd, jl, jk
1101      real zsqr, zalfa, zriw, zdel, zb, zalpha,zdz2n
1102      real zdelp, zdelpt
1103      REAL ZDZ2 (KLON,KLEV) , ZNORM(KLON) , zoro(KLON)
1104      REAL ZTAU (KLON,KLEV+1)
1105C
1106C-----------------------------------------------------------------------
1107C
1108C*         1.    INITIALIZATION
1109C                --------------
1110C
1111c      print *,' entree gwprofil'
1112 100  CONTINUE
1113C
1114C
1115C*    COMPUTATIONAL CONSTANTS.
1116C     ------------- ----------
1117C
1118      ilevh=KLEV/3
1119C
1120c     DO 400 ji=1,kgwd
1121c     jl=kdx(ji)
1122c  Modif vectorisation 02/04/2004
1123      DO 400 jl=kidia,kfdia
1124      if (ktest(jl).eq.1) then
1125      Zoro(JL)=Psig(JL)*Pdmod(JL)/4./max(pvar(jl),1.0)
1126      ZTAU(JL,KLEV+1)=PTAU(JL,KLEV+1)
1127      endif
1128  400 CONTINUE
1129 
1130C
1131      DO 430 JK=KLEV,2,-1
1132C
1133C
1134C*         4.1    CONSTANT WAVE STRESS UNTIL TOP OF THE
1135C                 BLOCKING LAYER.
1136  410 CONTINUE
1137C
1138c     DO 411 ji=1,kgwd
1139c     jl=kdx(ji)
1140c  Modif vectorisation 02/04/2004
1141      do 411 jl=kidia,kfdia
1142      if (ktest(jl).eq.1) then
1143           IF(JK.GT.KKCRITH(JL)) THEN
1144           PTAU(JL,JK)=ZTAU(JL,KLEV+1)
1145C          ENDIF
1146C          IF(JK.EQ.KKCRITH(JL)) THEN
1147           ELSE                   
1148           PTAU(JL,JK)=GRAHILO*ZTAU(JL,KLEV+1)
1149           ENDIF
1150      endif
1151 411  CONTINUE             
1152C
1153C*         4.15   CONSTANT SHEAR STRESS UNTIL THE TOP OF THE
1154C                 LOW LEVEL FLOW LAYER.
1155 415  CONTINUE
1156C       
1157C
1158C*         4.2    WAVE DISPLACEMENT AT NEXT LEVEL.
1159C
1160  420 CONTINUE
1161C
1162c     DO 421 ji=1,kgwd
1163c     jl=kdx(ji)
1164c  Modif vectorisation 02/04/2004
1165      do 421 jl=kidia,kfdia
1166      if(ktest(jl).eq.1) then
1167      IF(JK.LT.KKCRITH(JL)) THEN
1168      ZNORM(JL)=gkdrag*PRHO(JL,JK)*SQRT(PSTAB(JL,JK))*PVPH(JL,JK)
1169     *                                                    *zoro(jl)
1170      ZDZ2(JL,JK)=PTAU(JL,JK+1)/max(ZNORM(JL),gssec)
1171      ENDIF
1172      endif
1173  421 CONTINUE
1174C
1175C*         4.3    WAVE RICHARDSON NUMBER, NEW WAVE DISPLACEMENT
1176C*                AND STRESS:  BREAKING EVALUATION AND CRITICAL
1177C                 LEVEL
1178C
1179                         
1180c     DO 431 ji=1,kgwd
1181c     jl=Kdx(ji)
1182c  Modif vectorisation 02/04/2004
1183      do 431 jl=kidia,kfdia
1184      if(ktest(jl).eq.1) then
1185
1186          IF(JK.LT.KKCRITH(JL)) THEN
1187          IF((PTAU(JL,JK+1).LT.GTSEC).OR.(JK.LE.KCRIT(JL))) THEN
1188            PTAU(JL,JK)=0.0
1189          ELSE
1190               ZSQR=SQRT(PRI(JL,JK))
1191               ZALFA=SQRT(PSTAB(JL,JK)*ZDZ2(JL,JK))/PVPH(JL,JK)
1192               ZRIW=PRI(JL,JK)*(1.-ZALFA)/(1+ZALFA*ZSQR)**2
1193               IF(ZRIW.LT.GRCRIT) THEN
1194                 ZDEL=4./ZSQR/GRCRIT+1./GRCRIT**2+4./GRCRIT
1195                 ZB=1./GRCRIT+2./ZSQR
1196                 ZALPHA=0.5*(-ZB+SQRT(ZDEL))
1197                 ZDZ2N=(PVPH(JL,JK)*ZALPHA)**2/PSTAB(JL,JK)
1198                 PTAU(JL,JK)=ZNORM(JL)*ZDZ2N
1199               ELSE
1200                 PTAU(JL,JK)=ZNORM(JL)*ZDZ2(JL,JK)
1201               ENDIF
1202            PTAU(JL,JK)=MIN(PTAU(JL,JK),PTAU(JL,JK+1))
1203          ENDIF
1204          ENDIF
1205      endif
1206  431 CONTINUE
1207 
1208  430 CONTINUE
1209  440 CONTINUE
1210 
1211C  REORGANISATION OF THE STRESS PROFILE AT LOW LEVEL
1212
1213c     DO 530 ji=1,kgwd
1214c     jl=kdx(ji)
1215c  Modif vectorisation 02/04/2004
1216      do 530 jl=kidia,kfdia
1217      if(ktest(jl).eq.1) then
1218      ZTAU(JL,KKCRITH(JL))=PTAU(JL,KKCRITH(JL))
1219      ZTAU(JL,NSTRA)=PTAU(JL,NSTRA)
1220      endif
1221 530  CONTINUE     
1222
1223      DO 531 JK=1,KLEV
1224     
1225c     DO 532 ji=1,kgwd
1226c     jl=kdx(ji)
1227c  Modif vectorisation 02/04/2004
1228      do 532 jl=kidia,kfdia
1229      if(ktest(jl).eq.1) then
1230
1231               
1232         IF(JK.GT.KKCRITH(JL))THEN
1233
1234          ZDELP=PAPHM1(JL,JK)-PAPHM1(JL,KLEV+1    )
1235          ZDELPT=PAPHM1(JL,KKCRITH(JL))-PAPHM1(JL,KLEV+1    )
1236          PTAU(JL,JK)=ZTAU(JL,KLEV+1    ) +
1237     .                (ZTAU(JL,KKCRITH(JL))-ZTAU(JL,KLEV+1    ) )*
1238     .                ZDELP/ZDELPT
1239     
1240        ENDIF
1241           
1242      endif
1243 532  CONTINUE   
1244 
1245C  REORGANISATION IN THE STRATOSPHERE
1246
1247c     DO 533 ji=1,kgwd
1248c     jl=kdx(ji)
1249c  Modif vectorisation 02/04/2004
1250      do 533 jl=kidia,kfdia
1251      if(ktest(jl).eq.1) then
1252
1253
1254         IF(JK.LT.NSTRA)THEN
1255
1256          ZDELP =PAPHM1(JL,NSTRA)
1257          ZDELPT=PAPHM1(JL,JK)
1258          PTAU(JL,JK)=ZTAU(JL,NSTRA)*ZDELPT/ZDELP
1259
1260        ENDIF
1261
1262      endif
1263 533  CONTINUE
1264
1265C REORGANISATION IN THE TROPOSPHERE
1266
1267c      DO 534 ji=1,kgwd
1268c      jl=kdx(ji)
1269c  Modif vectorisation 02/04/2004
1270      do 534 jl=kidia,kfdia
1271      if(ktest(jl).eq.1) then
1272
1273
1274         IF(JK.LT.KKCRITH(JL).AND.JK.GT.NSTRA)THEN
1275
1276           ZDELP=PAPHM1(JL,JK)-PAPHM1(JL,KKCRITH(JL))
1277           ZDELPT=PAPHM1(JL,NSTRA)-PAPHM1(JL,KKCRITH(JL))
1278           PTAU(JL,JK)=ZTAU(JL,KKCRITH(JL)) +
1279     *                 (ZTAU(JL,NSTRA)-ZTAU(JL,KKCRITH(JL)))*ZDELP/ZDELPT
1280
1281       ENDIF
1282      endif
1283 534   CONTINUE
1284
1285 
1286 531  CONTINUE       
1287
1288
1289      RETURN
1290      END
1291      SUBROUTINE lift_noro (nlon,nlev,dtime,paprs,pplay,     
1292     e                   plat,pmea,pstd, ppic,
1293     e                   ktest,
1294     e                   t, u, v,
1295     s                   pulow, pvlow, pustr, pvstr,
1296     s                   d_t, d_u, d_v)
1297c
1298      IMPLICIT none
1299c======================================================================
1300c Auteur(s): F.Lott (LMD/CNRS) date: 19950201
1301c Objet: Frottement de la montagne Interface
1302c======================================================================
1303c Arguments:
1304c dtime---input-R- pas d'integration (s)
1305c paprs---input-R-pression pour chaque inter-couche (en Pa)
1306c pplay---input-R-pression pour le mileu de chaque couche (en Pa)
1307c t-------input-R-temperature (K)
1308c u-------input-R-vitesse horizontale (m/s)
1309c v-------input-R-vitesse horizontale (m/s)
1310c
1311c d_t-----output-R-increment de la temperature
1312c d_u-----output-R-increment de la vitesse u
1313c d_v-----output-R-increment de la vitesse v
1314c======================================================================
1315#include "dimensions.h"
1316#include "dimphy.h"
1317#include "YOMCST.h"
1318c
1319c ARGUMENTS
1320c
1321      INTEGER nlon,nlev
1322      REAL dtime
1323      REAL paprs(klon,klev+1)
1324      REAL pplay(klon,klev)
1325      REAL plat(nlon),pmea(nlon)
1326      REAL pstd(nlon)
1327      REAL ppic(nlon)
1328      REAL pulow(nlon),pvlow(nlon),pustr(nlon),pvstr(nlon)
1329      REAL t(nlon,nlev), u(nlon,nlev), v(nlon,nlev)
1330      REAL d_t(nlon,nlev), d_u(nlon,nlev), d_v(nlon,nlev)
1331c
1332      INTEGER i, k, ktest(nlon)
1333c
1334c Variables locales:
1335c
1336      REAL zgeom(klon,klev)
1337      REAL pdtdt(klon,klev), pdudt(klon,klev), pdvdt(klon,klev)
1338      REAL pt(klon,klev), pu(klon,klev), pv(klon,klev)
1339      REAL papmf(klon,klev),papmh(klon,klev+1)
1340c
1341c initialiser les variables de sortie (pour securite)
1342c
1343      DO i = 1,klon
1344         pulow(i) = 0.0
1345         pvlow(i) = 0.0
1346         pustr(i) = 0.0
1347         pvstr(i) = 0.0
1348      ENDDO
1349      DO k = 1, klev
1350      DO i = 1, klon
1351         d_t(i,k) = 0.0
1352         d_u(i,k) = 0.0
1353         d_v(i,k) = 0.0
1354         pdudt(i,k)=0.0
1355         pdvdt(i,k)=0.0
1356         pdtdt(i,k)=0.0
1357      ENDDO
1358      ENDDO
1359c
1360c preparer les variables d'entree (attention: l'ordre des niveaux
1361c verticaux augmente du haut vers le bas)
1362c
1363      DO k = 1, klev
1364      DO i = 1, klon
1365         pt(i,k) = t(i,klev-k+1)
1366         pu(i,k) = u(i,klev-k+1)
1367         pv(i,k) = v(i,klev-k+1)
1368         papmf(i,k) = pplay(i,klev-k+1)
1369      ENDDO
1370      ENDDO
1371      DO k = 1, klev+1
1372      DO i = 1, klon
1373         papmh(i,k) = paprs(i,klev-k+2)
1374      ENDDO
1375      ENDDO
1376      DO i = 1, klon
1377         zgeom(i,klev) = RD * pt(i,klev)
1378     .                  * LOG(papmh(i,klev+1)/papmf(i,klev))
1379      ENDDO
1380      DO k = klev-1, 1, -1
1381      DO i = 1, klon
1382         zgeom(i,k) = zgeom(i,k+1) + RD * (pt(i,k)+pt(i,k+1))/2.0
1383     .               * LOG(papmf(i,k+1)/papmf(i,k))
1384      ENDDO
1385      ENDDO
1386c
1387c appeler la routine principale
1388c
1389      CALL OROLIFT(klon,klev,ktest,
1390     .            dtime,
1391     .            papmh, zgeom,
1392     .            pt, pu, pv,
1393     .            plat,pmea, pstd, ppic,
1394     .            pulow,pvlow,
1395     .            pdudt,pdvdt,pdtdt)
1396C
1397      DO k = 1, klev
1398      DO i = 1, klon
1399         d_u(i,klev+1-k) = dtime*pdudt(i,k)
1400         d_v(i,klev+1-k) = dtime*pdvdt(i,k)
1401         d_t(i,klev+1-k) = dtime*pdtdt(i,k)
1402         pustr(i)        = pustr(i)
1403     .                    +RG*pdudt(i,k)*(papmh(i,k+1)-papmh(i,k))
1404         pvstr(i)        = pvstr(i)
1405     .                    +RG*pdvdt(i,k)*(papmh(i,k+1)-papmh(i,k))
1406      ENDDO
1407      ENDDO
1408c
1409      RETURN
1410      END
1411      SUBROUTINE OROLIFT( NLON,NLEV
1412     I                 , KTEST
1413     R                 , PTSPHY
1414     R                 , PAPHM1,PGEOM1,PTM1,PUM1,PVM1
1415     R                 , PLAT
1416     R                 , PMEA, PVAROR, ppic
1417C OUTPUTS
1418     R                 , PULOW,PVLOW
1419     R                 , PVOM,PVOL,PTE )
1420
1421C
1422C**** *OROLIFT: SIMULATE THE GEOSTROPHIC LIFT.
1423C
1424C     PURPOSE.
1425C     --------
1426C
1427C**   INTERFACE.
1428C     ----------
1429C          CALLED FROM *lift_noro
1430C     ----------
1431C
1432C     AUTHOR.
1433C     -------
1434C     F.LOTT  LMD 22/11/95
1435C
1436      implicit none
1437C
1438C
1439#include "dimensions.h"
1440#include "dimphy.h"
1441#include "YOMCST.h"
1442#include "YOEGWD.h"
1443C-----------------------------------------------------------------------
1444C
1445C*       0.1   ARGUMENTS
1446C              ---------
1447C
1448C
1449      integer nlon, nlev
1450      REAL  PTE(NLON,NLEV),
1451     *      PVOL(NLON,NLEV),
1452     *      PVOM(NLON,NLEV),
1453     *      PULOW(NLON),
1454     *      PVLOW(NLON)
1455      REAL  PUM1(NLON,NLEV),
1456     *      PVM1(NLON,NLEV),
1457     *      PTM1(NLON,NLEV),
1458     *      PLAT(NLON),PMEA(NLON),
1459     *      PVAROR(NLON),
1460     *      ppic(NLON),
1461     *      PGEOM1(NLON,NLEV),
1462     *      PAPHM1(NLON,NLEV+1)
1463C
1464      INTEGER  KTEST(NLON)
1465      real ptsphy
1466C-----------------------------------------------------------------------
1467C
1468C*       0.2   LOCAL ARRAYS
1469C              ------------
1470      logical lifthigh, ll1
1471      integer klevm1, jl, ilevh, jk
1472      real zcons1, ztmst, zrtmst,zpi, zhgeo
1473      real zdelp, zslow, zsqua, zscav, zbet
1474      INTEGER 
1475     *         IKNUB(klon),
1476     *         IKNUL(klon)
1477      LOGICAL LL1(KLON,KLEV+1)
1478C
1479      REAL   ZTAU(KLON,KLEV+1),
1480     *       ZTAV(KLON,KLEV+1),
1481     *       ZRHO(KLON,KLEV+1)
1482      REAL   ZDUDT(KLON),
1483     *       ZDVDT(KLON)
1484      REAL ZHCRIT(KLON,KLEV)
1485C-----------------------------------------------------------------------
1486C
1487C*         1.1  INITIALIZATIONS
1488C               ---------------
1489
1490      LIFTHIGH=.FALSE.
1491
1492      IF(NLON.NE.KLON.OR.NLEV.NE.KLEV)STOP
1493      ZCONS1=1./RD
1494      KLEVM1=KLEV-1
1495      ZTMST=PTSPHY
1496      ZRTMST=1./ZTMST
1497      ZPI=ACOS(-1.)
1498C
1499      DO 1001 JL=kidia,kfdia
1500      ZRHO(JL,KLEV+1)  =0.0
1501      PULOW(JL)        =0.0
1502      PVLOW(JL)        =0.0
1503      iknub(JL)   =klev
1504      iknul(JL)   =klev
1505      ilevh=klev/3
1506      ll1(jl,klev+1)=.false.
1507      DO 1000 JK=1,KLEV
1508      PVOM(JL,JK)=0.0
1509      PVOL(JL,JK)=0.0
1510      PTE (JL,JK)=0.0
1511 1000 CONTINUE
1512 1001 CONTINUE
1513
1514C
1515C*         2.1     DEFINE LOW LEVEL WIND, PROJECT WINDS IN PLANE OF
1516C*                 LOW LEVEL WIND, DETERMINE SECTOR IN WHICH TO TAKE
1517C*                 THE VARIANCE AND SET INDICATOR FOR CRITICAL LEVELS.
1518C
1519C
1520C
1521      DO 2006 JK=KLEV,1,-1
1522      DO 2007 JL=kidia,kfdia
1523      IF(KTEST(JL).EQ.1) THEN
1524      ZHCRIT(JL,JK)=amax1(Ppic(JL)-pmea(JL),100.)
1525      ZHGEO=PGEOM1(JL,JK)/RG
1526      ll1(JL,JK)=(ZHGEO.GT.ZHCRIT(JL,JK))
1527      IF(ll1(JL,JK).XOR.ll1(JL,JK+1)) THEN
1528        iknub(JL)=JK
1529      ENDIF
1530      ENDIF
1531 2007 CONTINUE
1532 2006 CONTINUE
1533C
1534      do 2010 jl=kidia,kfdia
1535      IF(KTEST(JL).EQ.1) THEN
1536      iknub(jl)=max(iknub(jl),klev/2)
1537      iknul(jl)=max(iknul(jl),2*klev/3)
1538      if(iknub(jl).gt.nktopg) iknub(jl)=nktopg
1539      if(iknub(jl).eq.nktopg) iknul(jl)=klev
1540      if(iknub(jl).eq.iknul(jl)) iknub(jl)=iknul(jl)-1
1541      ENDIF
1542 2010 continue
1543
1544C     do 2011 jl=kidia,kfdia
1545C     IF(KTEST(JL).EQ.1) THEN
1546C       print *,' iknul= ',iknul(jl),'  iknub=',iknub(jl)
1547C     ENDIF
1548C2011 continue
1549
1550C     PRINT *,'  DANS OROLIFT: 2010'
1551
1552      DO 223 JK=KLEV,2,-1
1553      DO 222 JL=kidia,kfdia
1554        ZRHO(JL,JK)=2.*PAPHM1(JL,JK)*ZCONS1/(PTM1(JL,JK)+PTM1(JL,JK-1))
1555  222 CONTINUE
1556  223 CONTINUE
1557C     PRINT *,'  DANS OROLIFT: 223'
1558
1559C********************************************************************
1560C
1561C*     DEFINE LOW LEVEL FLOW
1562C      -------------------
1563      DO 2115 JK=klev,1,-1
1564      DO 2116 JL=kidia,kfdia
1565      IF(KTEST(JL).EQ.1) THEN
1566      if(jk.ge.iknub(jl).and.jk.le.iknul(jl)) then
1567        pulow(JL)=pulow(JL)+PUM1(JL,JK)*(PAPHM1(JL,JK+1)-PAPHM1(JL,JK))
1568        pvlow(JL)=pvlow(JL)+PVM1(JL,JK)*(PAPHM1(JL,JK+1)-PAPHM1(JL,JK))
1569        zrho(JL,klev+1)=zrho(JL,klev+1)
1570     *                 +zrho(JL,JK)*(PAPHM1(JL,JK+1)-PAPHM1(JL,JK))
1571      end if
1572      ENDIF
1573 2116 CONTINUE
1574 2115 CONTINUE
1575      DO 2110 JL=kidia,kfdia
1576      IF(KTEST(JL).EQ.1) THEN
1577      pulow(JL)=pulow(JL)/(PAPHM1(JL,iknul(jl)+1)-PAPHM1(JL,iknub(jl)))
1578      pvlow(JL)=pvlow(JL)/(PAPHM1(JL,iknul(jl)+1)-PAPHM1(JL,iknub(jl)))
1579      zrho(JL,klev+1)=zrho(JL,klev+1)
1580     *               /(PAPHM1(JL,iknul(jl)+1)-PAPHM1(JL,iknub(jl)))
1581      ENDIF
1582 2110 CONTINUE
1583
1584
1585200   CONTINUE
1586
1587C***********************************************************
1588C
1589C*         3.      COMPUTE MOUNTAIN LIFT
1590C
1591  300 CONTINUE
1592C
1593      DO 301 JL=kidia,kfdia
1594      IF(KTEST(JL).EQ.1) THEN
1595       ZTAU(JL,KLEV+1)= - GKLIFT*ZRHO(JL,KLEV+1)*2.*ROMEGA*
1596C    *                 (2*PVAROR(JL)+PMEA(JL))*
1597     *                 2*PVAROR(JL)*
1598     *                 SIN(ZPI/180.*PLAT(JL))*PVLOW(JL)
1599       ZTAV(JL,KLEV+1)=   GKLIFT*ZRHO(JL,KLEV+1)*2.*ROMEGA*
1600C    *                 (2*PVAROR(JL)+PMEA(JL))*
1601     *                 2*PVAROR(JL)*
1602     *                 SIN(ZPI/180.*PLAT(JL))*PULOW(JL)
1603      ELSE
1604       ZTAU(JL,KLEV+1)=0.0
1605       ZTAV(JL,KLEV+1)=0.0
1606      ENDIF
1607301   CONTINUE
1608
1609C
1610C*         4.      COMPUTE LIFT PROFILE         
1611C*                 --------------------   
1612C
1613
1614  400 CONTINUE
1615
1616      DO 401 JK=1,KLEV
1617      DO 401 JL=kidia,kfdia
1618      IF(KTEST(JL).EQ.1) THEN
1619      ZTAU(JL,JK)=ZTAU(JL,KLEV+1)*PAPHM1(JL,JK)/PAPHM1(JL,KLEV+1)
1620      ZTAV(JL,JK)=ZTAV(JL,KLEV+1)*PAPHM1(JL,JK)/PAPHM1(JL,KLEV+1)
1621      ELSE
1622      ZTAU(JL,JK)=0.0
1623      ZTAV(JL,JK)=0.0
1624      ENDIF
1625401   CONTINUE
1626C
1627C
1628C*         5.      COMPUTE TENDENCIES.
1629C*                 -------------------
1630      IF(LIFTHIGH)THEN
1631C
1632  500 CONTINUE
1633C     PRINT *,'  DANS OROLIFT: 500'
1634C
1635C  EXPLICIT SOLUTION AT ALL LEVELS
1636C
1637      DO 524 JK=1,klev
1638      DO 523 JL=KIDIA,KFDIA
1639      IF(KTEST(JL).EQ.1) THEN
1640      ZDELP=PAPHM1(JL,JK+1)-PAPHM1(JL,JK)
1641      ZDUDT(JL)=-RG*(ZTAU(JL,JK+1)-ZTAU(JL,JK))/ZDELP
1642      ZDVDT(JL)=-RG*(ZTAV(JL,JK+1)-ZTAV(JL,JK))/ZDELP
1643      ENDIF 
1644  523 CONTINUE
1645  524 CONTINUE
1646C
1647C  PROJECT PERPENDICULARLY TO U NOT TO DESTROY ENERGY
1648C
1649      DO 530 JK=1,klev
1650      DO 530 JL=KIDIA,KFDIA
1651      IF(KTEST(JL).EQ.1) THEN
1652
1653        ZSLOW=SQRT(PULOW(JL)**2+PVLOW(JL)**2)
1654        ZSQUA=AMAX1(SQRT(PUM1(JL,JK)**2+PVM1(JL,JK)**2),GVSEC)
1655        ZSCAV=-ZDUDT(JL)*PVM1(JL,JK)+ZDVDT(JL)*PUM1(JL,JK)
1656        IF(ZSQUA.GT.GVSEC)THEN
1657          PVOM(JL,JK)=-ZSCAV*PVM1(JL,JK)/ZSQUA**2
1658          PVOL(JL,JK)= ZSCAV*PUM1(JL,JK)/ZSQUA**2
1659        ELSE
1660          PVOM(JL,JK)=0.0
1661          PVOL(JL,JK)=0.0     
1662        ENDIF 
1663        ZSQUA=SQRT(PUM1(JL,JK)**2+PUM1(JL,JK)**2)               
1664        IF(ZSQUA.LT.ZSLOW)THEN
1665          PVOM(JL,JK)=ZSQUA/ZSLOW*PVOM(JL,JK)
1666          PVOL(JL,JK)=ZSQUA/ZSLOW*PVOL(JL,JK)
1667        ENDIF
1668
1669      ENDIF 
1670530   CONTINUE
1671C
1672C  6.  LOW LEVEL LIFT, SEMI IMPLICIT:
1673C  ----------------------------------
1674
1675      ELSE
1676
1677        DO 601 JL=KIDIA,KFDIA
1678        IF(KTEST(JL).EQ.1) THEN
1679          DO JK=KLEV,IKNUB(JL),-1
1680          ZBET=GKLIFT*2.*ROMEGA*SIN(ZPI/180.*PLAT(JL))*ztmst*
1681     *        (PGEOM1(JL,IKNUB(JL)-1)-PGEOM1(JL,  JK))/
1682     *        (PGEOM1(JL,IKNUB(JL)-1)-PGEOM1(JL,KLEV))
1683          ZDUDT(JL)=-PUM1(JL,JK)/ztmst/(1+ZBET**2)
1684          ZDVDT(JL)=-PVM1(JL,JK)/ztmst/(1+ZBET**2)
1685          PVOM(JL,JK)= ZBET**2*ZDUDT(JL) - ZBET   *ZDVDT(JL)
1686          PVOL(JL,JK)= ZBET   *ZDUDT(JL) + ZBET**2*ZDVDT(JL)   
1687          ENDDO
1688        ENDIF
1689 601    CONTINUE
1690
1691      ENDIF
1692
1693      RETURN
1694      END
1695      SUBROUTINE SUGWD(NLON,NLEV,paprs,pplay)
1696C
1697C**** *SUGWD* INITIALIZE COMMON YOEGWD CONTROLLING GRAVITY WAVE DRAG
1698C
1699C     PURPOSE.
1700C     --------
1701C           INITIALIZE YOEGWD, THE COMMON THAT CONTROLS THE
1702C           GRAVITY WAVE DRAG PARAMETRIZATION.
1703C
1704C**   INTERFACE.
1705C     ----------
1706C        CALL *SUGWD* FROM *SUPHEC*
1707C              -----        ------
1708C
1709C        EXPLICIT ARGUMENTS :
1710C        --------------------
1711C        PSIG        : VERTICAL COORDINATE TABLE
1712C        NLEV        : NUMBER OF MODEL LEVELS
1713C
1714C        IMPLICIT ARGUMENTS :
1715C        --------------------
1716C        COMMON YOEGWD
1717C
1718C     METHOD.
1719C     -------
1720C        SEE DOCUMENTATION
1721C
1722C     EXTERNALS.
1723C     ----------
1724C        NONE
1725C
1726C     REFERENCE.
1727C     ----------
1728C        ECMWF Research Department documentation of the IFS
1729C
1730C     AUTHOR.
1731C     -------
1732C        MARTIN MILLER             *ECMWF*
1733C
1734C     MODIFICATIONS.
1735C     --------------
1736C        ORIGINAL : 90-01-01
1737C     ------------------------------------------------------------------
1738      implicit none
1739C
1740C     -----------------------------------------------------------------
1741#include "YOEGWD.h"
1742C      ----------------------------------------------------------------
1743C
1744      integer nlon,nlev, jk
1745      REAL paprs(nlon,nlev+1)
1746      REAL pplay(nlon,nlev)
1747      real zpr,zstra,zsigt,zpm1r
1748C
1749C*       1.    SET THE VALUES OF THE PARAMETERS
1750C              --------------------------------
1751C
1752 100  CONTINUE
1753C
1754      PRINT *,' DANS SUGWD NLEV=',NLEV
1755      GHMAX=10000.
1756C
1757      ZPR=100000.
1758      ZSTRA=0.1
1759      ZSIGT=0.94
1760cold  ZPR=80000.
1761cold  ZSIGT=0.85
1762C
1763      DO 110 JK=1,NLEV
1764      ZPM1R=pplay(nlon/2,jk)/paprs(nlon/2,1)
1765      IF(ZPM1R.GE.ZSIGT)THEN
1766         nktopg=JK
1767      ENDIF
1768      ZPM1R=pplay(nlon/2,jk)/paprs(nlon/2,1)
1769      IF(ZPM1R.GE.ZSTRA)THEN
1770         NSTRA=JK
1771      ENDIF
1772  110 CONTINUE
1773c
1774c  inversion car dans orodrag on compte les niveaux a l'envers
1775      nktopg=nlev-nktopg+1
1776      nstra=nlev-nstra
1777      print *,' DANS SUGWD nktopg=', nktopg
1778      print *,' DANS SUGWD nstra=', nstra
1779C
1780      GSIGCR=0.80
1781C
1782      GKDRAG=0.2
1783      GRAHILO=1.   
1784      GRCRIT=0.01
1785      GFRCRIT=1.0
1786      GKWAKE=0.50
1787C
1788      GKLIFT=0.50 
1789      GVCRIT =0.0
1790C
1791C
1792C      ----------------------------------------------------------------
1793C
1794C*       2.    SET VALUES OF SECURITY PARAMETERS
1795C              ---------------------------------
1796C
1797 200  CONTINUE
1798C
1799      GVSEC=0.10
1800      GSSEC=1.E-12
1801C
1802      GTSEC=1.E-07
1803C
1804C      ----------------------------------------------------------------
1805C
1806      RETURN
1807      END
Note: See TracBrowser for help on using the repository browser.