source: LMDZ5/branches/IPSLCM6.0.11pre/libf/phylmd/stratosphere_mask.F90 @ 5360

Last change on this file since 5360 was 2773, checked in by Laurent Fairhead, 8 years ago

ooops
LF

  • Property svn:keywords set to Id
File size: 7.5 KB
Line 
1!
2! $Id: stratosphere_mask.F90 2773 2017-01-25 08:41:17Z abarral $
3!
4SUBROUTINE stratosphere_mask(t_seri, pplay, xlat)
5
6!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
7!
8! determination of tropopause height and temperature from gridded temperature data
9!
10! reference: Reichler, T., M. Dameris, and R. Sausen (2003):
11! modified: 6/28/06 tjr
12! adapted to LMDZ by C. Kleinschmitt (2016-02-15)
13!
14! input:    temp(nlon,nlat,nlev)    3D-temperature field
15!       ps(nlon,nlat)       2D-surface pressure field
16!       zs(nlon,nlat)       2D-surface height
17!       nlon            grid points in x
18!       nlat            grid points in y
19!       pfull(nlon,nlat,nlev)   full pressure levels in Pa
20!       plimu           upper limit for tropopause pressure
21!       pliml           lower limit for tropopause pressure
22!       gamma           tropopause criterion, e.g. -0.002 K/m
23!
24! output:   tp(nlon, nlat)      tropopause pressure in Pa, -999. if undefined
25!           ttp(nlon, nlat)     tropopause temperature in K, -999. if undefined
26!           ztp(nlon, nlat)     tropopause height in m, -999. if undefined
27!       tperr           # of undetermined values
28!
29!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
30
31USE dimphy
32USE phys_local_var_mod, ONLY: stratomask
33#ifdef CPP_StratAer
34USE phys_local_var_mod, ONLY: p_tropopause
35#endif
36USE print_control_mod, ONLY: lunout, prt_level
37
38IMPLICIT NONE
39
40REAL,DIMENSION(klon,klev),INTENT(IN)   :: t_seri  ! Temperature
41REAL,DIMENSION(klon,klev),INTENT(IN)   :: pplay   ! pression pour le mileu de chaque couche (en Pa)
42REAL,DIMENSION(klon),INTENT(IN)        :: xlat    ! latitudes pour chaque point
43
44REAL, PARAMETER                        :: plimu=45000.
45REAL, PARAMETER                        :: pliml=7500.
46REAL, PARAMETER                        :: gamma=-0.002
47LOGICAL, PARAMETER                     :: dofill=.true.
48REAL,DIMENSION(klon)                   :: tp, ttp, ztp
49REAL,DIMENSION(klev)                   :: t, p
50INTEGER                                :: tperr, i, k, invert, ifil
51REAL                                   :: ptrp, ttrp, ztrp, psrf, zsrf, pi
52
53pi     = 4.*ATAN(1.)
54
55tperr = 0
56DO i=1,klon
57  DO k=1,klev
58    t(k)=t_seri(i,klev+1-k)
59    p(k)=pplay(i,klev+1-k)
60  ENDDO
61  psrf=pplay(i,1)
62  zsrf=0.0
63  call twmo(klev, t, p, psrf, zsrf, plimu, pliml, gamma, ptrp, ttrp, ztrp)
64  tp(i)=ptrp
65  ttp(i)=ttrp
66  ztp(i)=ztrp
67  IF (ptrp.lt.0.0) THEN
68    tperr = tperr+1
69  ENDIF
70ENDDO
71
72! fill holes
73IF (dofill) THEN
74  ifil=0
75  DO i=1,klon
76  IF (tp(i).lt.-990.) THEN
77    !set missing values to very simple profile (neighbour averaging too expensive in LMDZ)
78    tp(i)=50000.-20000.*cos(xlat(i)/360.*2.*pi)
79    ifil=ifil+1
80  ENDIF
81  ENDDO
82!
83  IF (ifil.ne.tperr) THEN
84    CALL abort_physic('stratosphere_mask', 'inconsistency',1)
85  ENDIF
86ENDIF
87!
88DO i=1, klon
89DO k=1, klev
90  IF (pplay(i,k).LT.tp(i)) THEN
91    stratomask(i,k)=1.0
92  ELSE
93    stratomask(i,k)=0.0
94  ENDIF
95ENDDO
96ENDDO
97
98!--this is only diagnosedd in the case of StratAer
99!--but it could be useful to LMDz
100#ifdef CPP_StratAer
101p_tropopause(:)=tp(:)
102#endif
103
104IF (ifil.gt.0 .and. prt_level >5) THEN
105  write(lunout,*)'Tropopause: number of undetermined values =', ifil
106ENDIF
107
108RETURN
109END SUBROUTINE stratosphere_mask
110
111!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
112! twmo
113!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
114
115subroutine twmo(level, t, p, ps, zs, plimu, pliml, gamma, ptrp, ttrp, ztrp)
116
117implicit none
118
119include "YOMCST.h"
120
121integer,intent(in)              :: level
122real,intent(in),dimension(level):: t, p
123real,intent(in)                 :: plimu, pliml, gamma, ps, zs
124real,intent(out)                :: ptrp, ttrp, ztrp
125
126real,parameter                  :: deltaz = 2000.0
127
128real     :: faktor
129real     :: pmk, pm, a, b, tm, dtdp, dtdz, dlnp, tdlnp
130real     :: ag, bg, ptph, ttph, a0, b0
131real     :: pm0, tm0, pmk0, dtdz0
132real     :: p2km, asum, aquer
133real     :: pmk2, pm2, a2, b2, tm2, dtdp2, dtdz2
134integer  :: icount, jj, j
135
136ptrp=-999.0             
137ttrp=-999.0                 
138ztrp=-999.0                 
139
140faktor = -RG/R
141
142do j=level,2,-1
143
144   ! dt/dz
145   pmk= 0.5 * (p(j-1)**rkappa+p(j)**rkappa)        ! p**k at half level
146   pm = pmk**(1./rkappa)                    ! p at half level
147   a = (t(j-1)-t(j))/(p(j-1)**rkappa-p(j)**rkappa)    ! dT/dp^k
148   b = t(j)-(a*p(j)**rkappa)
149   tm = a * pmk + b                    ! T at half level         
150   dtdp = a * rkappa * pm**(rkappa-1.) ! dT/dp at half level
151   dtdz = faktor*dtdp*pm/tm            ! dT/dz at half level
152
153   ! dt/dz valid?
154   if (j.eq.level)     go to 999                 ! no, start level, initialize first
155   if (dtdz.le.gamma)  go to 999                 ! no, dt/dz < -2 K/km
156   if (dtdz0.gt.gamma) go to 999                 ! no, dt/dz below > -2 K/km
157   if (pm.gt.plimu)    go to 999                 ! no, pm too low
158
159   ! dtdz is valid, calculate tropopause pressure ptph
160   ag = (dtdz-dtdz0) / (pmk-pmk0)     
161   bg = dtdz0 - (ag * pmk0)         
162   ptph = exp(log((gamma-bg)/ag)/rkappa)
163
164   ! calculate temperature at this ptph assuming linear gamma
165   ! interpolation
166   ttph = tm0
167   ttph = ttph - (bg * log(pm0)  + ag * (pm0**rkappa) /rkappa) / faktor*t(j)
168   ttph = ttph + (bg * log(ptph) + ag * (ptph**rkappa)/rkappa) / faktor*t(j)
169
170   if (ptph.lt.pliml) go to 999
171   if (ptph.gt.plimu) go to 999     
172
173   ! 2nd test: dtdz above 2 km must not exceed gamma
174   p2km = ptph + deltaz*(pm/tm)*faktor      ! p at ptph + 2km
175   asum = 0.0                               ! dtdz above
176   icount = 0                               ! number of levels above
177
178   ! test until apm < p2km
179   do jj=j,2,-1
180
181       pmk2 = .5 * (p(jj-1)**rkappa+p(jj)**rkappa)    ! p mean ^kappa
182       pm2 = pmk2**(1/rkappa)                      ! p mean
183       if(pm2.gt.ptph) go to 110                ! doesn't happen
184       if(pm2.lt.p2km) go to 888                ! ptropo is valid
185
186       a2 = (t(jj-1)-t(jj))                     ! a
187       a2 = a2/(p(jj-1)**rkappa-p(jj)**rkappa)
188       b2 = t(jj)-(a2*p(jj)**rkappa)               ! b
189       tm2 = a2 * pmk2 + b2                     ! T mean
190       dtdp2 = a2 * rkappa * (pm2**(rkappa-1))        ! dt/dp
191       dtdz2 = faktor*dtdp2*pm2/tm2
192       asum = asum+dtdz2
193       icount = icount+1
194       aquer = asum/float(icount)               ! dt/dz mean
195   
196       ! discard ptropo ?
197        if (aquer.le.gamma) go to 999           ! dt/dz above < gamma
198
199110 continue
200    enddo                   ! test next level
201
202888 continue                    ! ptph is valid
203    ptrp = ptph
204    ttrp = ttph
205
206! now calculate height of tropopause by integrating hypsometric equation
207! from ps to ptrp
208! linearly interpolate in p (results are identical to p^kappa)
209
210jj = LEVEL                                       ! bottom
211do while ((P(jj).gt.PS) .or. (T(jj).lt.100))     ! T must be valid too
212  jj=jj-1
213enddo
214
215DLNP = log(PS/P(jj))                             ! from surface pressure
216TM = T(jj)                                       ! take TM of lowest level (better: extrapolate)
217TDLNP = TM*DLNP
218
219do while ( (JJ.ge.2) .and. (PTRP.lt.P(jj-1)) )
220  DLNP = log(P(jj)/P(jj-1))
221  TM = 0.5 * (T(jj) + T(jj-1))
222  TDLNP = TDLNP + TM*DLNP
223  JJ=JJ-1
224enddo
225
226DLNP = log(P(jj)/PTRP)                           ! up to tropopause pressure
227TM = 0.5 * (T(jj) + TTRP)                        ! use TTRP to get TM of this level
228TDLNP = TDLNP + TM*DLNP
229
230ZTRP = ZS + TDLNP*R/RG
231
232if (ZTRP .lt. 0) then
233  print*,ZTRP
234  print*,PS
235  print*,P
236  print*,T
237  print*,ZS
238  stop
239endif
240
241return
242
243999 continue                    ! continue search at next higher level
244    tm0 = tm
245    pm0 = pm
246    pmk0 = pmk
247    dtdz0  = dtdz
248    a0 = a
249    b0 = b
250
251enddo
252
253! no tropopouse found
254return
255end subroutine twmo
Note: See TracBrowser for help on using the repository browser.