source: LMDZ6/branches/Ocean_skin/libf/phylmd/stratosphere_mask.F90 @ 5434

Last change on this file since 5434 was 3123, checked in by oboucher, 7 years ago

Compute z_tropopause for diagnostic purpose from geoid rather than from the surface as per CF convention

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