1 | ! |
---|
2 | ! $Id: stratosphere_mask.F90 2773 2017-01-25 08:41:17Z aclsce $ |
---|
3 | ! |
---|
4 | SUBROUTINE 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 | |
---|
31 | USE dimphy |
---|
32 | USE phys_local_var_mod, ONLY: stratomask |
---|
33 | #ifdef CPP_StratAer |
---|
34 | USE phys_local_var_mod, ONLY: p_tropopause |
---|
35 | #endif |
---|
36 | USE print_control_mod, ONLY: lunout, prt_level |
---|
37 | |
---|
38 | IMPLICIT NONE |
---|
39 | |
---|
40 | REAL,DIMENSION(klon,klev),INTENT(IN) :: t_seri ! Temperature |
---|
41 | REAL,DIMENSION(klon,klev),INTENT(IN) :: pplay ! pression pour le mileu de chaque couche (en Pa) |
---|
42 | REAL,DIMENSION(klon),INTENT(IN) :: xlat ! latitudes pour chaque point |
---|
43 | |
---|
44 | REAL, PARAMETER :: plimu=45000. |
---|
45 | REAL, PARAMETER :: pliml=7500. |
---|
46 | REAL, PARAMETER :: gamma=-0.002 |
---|
47 | LOGICAL, PARAMETER :: dofill=.true. |
---|
48 | REAL,DIMENSION(klon) :: tp, ttp, ztp |
---|
49 | REAL,DIMENSION(klev) :: t, p |
---|
50 | INTEGER :: tperr, i, k, invert, ifil |
---|
51 | REAL :: ptrp, ttrp, ztrp, psrf, zsrf, pi |
---|
52 | |
---|
53 | pi = 4.*ATAN(1.) |
---|
54 | |
---|
55 | tperr = 0 |
---|
56 | DO 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 |
---|
70 | ENDDO |
---|
71 | |
---|
72 | ! fill holes |
---|
73 | IF (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 |
---|
86 | ENDIF |
---|
87 | ! |
---|
88 | DO i=1, klon |
---|
89 | DO 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 |
---|
95 | ENDDO |
---|
96 | ENDDO |
---|
97 | |
---|
98 | !--this is only diagnosedd in the case of StratAer |
---|
99 | !--but it could be useful to LMDz |
---|
100 | #ifdef CPP_StratAer |
---|
101 | p_tropopause(:)=tp(:) |
---|
102 | #endif |
---|
103 | |
---|
104 | IF (ifil.gt.0 .and. prt_level >5) THEN |
---|
105 | write(lunout,*)'Tropopause: number of undetermined values =', ifil |
---|
106 | ENDIF |
---|
107 | |
---|
108 | RETURN |
---|
109 | END SUBROUTINE stratosphere_mask |
---|
110 | |
---|
111 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
112 | ! twmo |
---|
113 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
114 | |
---|
115 | subroutine twmo(level, t, p, ps, zs, plimu, pliml, gamma, ptrp, ttrp, ztrp) |
---|
116 | |
---|
117 | implicit none |
---|
118 | |
---|
119 | include "YOMCST.h" |
---|
120 | |
---|
121 | integer,intent(in) :: level |
---|
122 | real,intent(in),dimension(level):: t, p |
---|
123 | real,intent(in) :: plimu, pliml, gamma, ps, zs |
---|
124 | real,intent(out) :: ptrp, ttrp, ztrp |
---|
125 | |
---|
126 | real,parameter :: deltaz = 2000.0 |
---|
127 | |
---|
128 | real :: faktor |
---|
129 | real :: pmk, pm, a, b, tm, dtdp, dtdz, dlnp, tdlnp |
---|
130 | real :: ag, bg, ptph, ttph, a0, b0 |
---|
131 | real :: pm0, tm0, pmk0, dtdz0 |
---|
132 | real :: p2km, asum, aquer |
---|
133 | real :: pmk2, pm2, a2, b2, tm2, dtdp2, dtdz2 |
---|
134 | integer :: icount, jj, j |
---|
135 | |
---|
136 | ptrp=-999.0 |
---|
137 | ttrp=-999.0 |
---|
138 | ztrp=-999.0 |
---|
139 | |
---|
140 | faktor = -RG/R |
---|
141 | |
---|
142 | do 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 | |
---|
199 | 110 continue |
---|
200 | enddo ! test next level |
---|
201 | |
---|
202 | 888 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 | |
---|
210 | jj = LEVEL ! bottom |
---|
211 | do while ((P(jj).gt.PS) .or. (T(jj).lt.100)) ! T must be valid too |
---|
212 | jj=jj-1 |
---|
213 | enddo |
---|
214 | |
---|
215 | DLNP = log(PS/P(jj)) ! from surface pressure |
---|
216 | TM = T(jj) ! take TM of lowest level (better: extrapolate) |
---|
217 | TDLNP = TM*DLNP |
---|
218 | |
---|
219 | do 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 |
---|
224 | enddo |
---|
225 | |
---|
226 | DLNP = log(P(jj)/PTRP) ! up to tropopause pressure |
---|
227 | TM = 0.5 * (T(jj) + TTRP) ! use TTRP to get TM of this level |
---|
228 | TDLNP = TDLNP + TM*DLNP |
---|
229 | |
---|
230 | ZTRP = ZS + TDLNP*R/RG |
---|
231 | |
---|
232 | if (ZTRP .lt. 0) then |
---|
233 | print*,ZTRP |
---|
234 | print*,PS |
---|
235 | print*,P |
---|
236 | print*,T |
---|
237 | print*,ZS |
---|
238 | stop |
---|
239 | endif |
---|
240 | |
---|
241 | return |
---|
242 | |
---|
243 | 999 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 | |
---|
251 | enddo |
---|
252 | |
---|
253 | ! no tropopouse found |
---|
254 | return |
---|
255 | end subroutine twmo |
---|