3 |
michaesp |
1 |
PROGRAM clustering
|
|
|
2 |
|
|
|
3 |
! *********************************************************************************
|
|
|
4 |
! * Interpolate fields to tropopause and assign labels with clustering algorithm *
|
|
|
5 |
! * Michael Sprenger / Spring 2006: initial version *
|
|
|
6 |
! * Michael Sprenger / 2012: horizontal propagation of LABEL 2 prohibited in *
|
|
|
7 |
! * lower troposphere *
|
|
|
8 |
! * Bojan Skerlak / Spring 2012 / Fall 2012: advanced LABEL 1-5 handling *
|
|
|
9 |
! * Jan 2013: - change horiz. prop. decision method to height(log(pressure)) *
|
|
|
10 |
! * - perform check of label for exotic cases where the pressure *
|
|
|
11 |
! * is higher than 800 hPa (deep folds) *
|
|
|
12 |
! * - additional cutoff-special-treatment prohibits identification *
|
|
|
13 |
! * of strat. air as cutoffs if below horizontal prop. limit *
|
|
|
14 |
! * (e.g. deep 'funnels' with 'bumps' on the side) *
|
|
|
15 |
! * - additional check in strat. (TH>380) for 'wrong' sign of PV *
|
|
|
16 |
! * that can lead to errors in calculating TROPO *
|
|
|
17 |
! * - additional check for PV towers (whole column gets label 2) *
|
|
|
18 |
! * July 2013: - replace fixed tropo_pv,tropo_th by user-defined values *
|
|
|
19 |
! * Jan 2014: - v2.0 Add specific humidity criterion for diabatic PV *
|
|
|
20 |
! * Feb 2014: - v2.1 Correct 'filling' when labels 2/3 or 2/5 meet *
|
|
|
21 |
! * Mar 2014: - v2.2 Filling of 2/3 and 2/5 also propagates downward via *
|
|
|
22 |
! * 'connect' subroutine. Thus, now also dird/diru have to be *
|
|
|
23 |
! * specified. *
|
|
|
24 |
! * - v2.3 Identify and write out tropopause folds (FOLD) 1/0 *
|
|
|
25 |
! * Nov 2014: - remove some commented sections and delete temporary stuff *
|
|
|
26 |
! * Dec 2014: - simply the code to only contain the core algorithm and no *
|
|
|
27 |
! * input/output section (which depends on data format and *
|
|
|
28 |
! * model setup) *
|
|
|
29 |
! * Andrea Pozzer (MPIC):
|
|
|
30 |
! * Dec 2014: - Implementation for netcdf I/O, Makefile and namelist reading *
|
|
|
31 |
! *********************************************************************************
|
|
|
32 |
|
|
|
33 |
!????????????????????????????????????????????????????????????????????????????????
|
|
|
34 |
!????????????????????????????????????????????????????????????????????????????????
|
|
|
35 |
!????????????????????????????????????????????????????????????????????????????????
|
|
|
36 |
!??????????????? ??????????????????????????
|
|
|
37 |
!??????????????? 1) 3D-LABELLING ??????????????????????????
|
|
|
38 |
!??????????????? 2) TROPOPAUSE FOLD IDENTIFICATION ??????????????????????????
|
|
|
39 |
!??????????????? ??????????????????????????
|
|
|
40 |
!????????????????????????????????????????????????????????????????????????????????
|
|
|
41 |
!????????????????????????????????????????????????????????????????????????????????
|
|
|
42 |
!MMM?????????????????????????????????????????????????????????????????????????????
|
|
|
43 |
!.....DMMI???????????????????????????????????????????MMMMMMMMMMM?I???????????????
|
|
|
44 |
!.........IMO????????????????????????????????????IMM ........... 8MM8???????????
|
|
|
45 |
!............MM?????????????????????????????????M$ ..............MMMI??????
|
|
|
46 |
!............. MD?????????????????????????????IM. PMIN .................OMMI??
|
|
|
47 |
!................MI???????????????????????????M . ^ ......................NM
|
|
|
48 |
!.................MO?????????????????????????M... | .............................
|
|
|
49 |
!.................. M???????????????????????M.... | .............................
|
|
|
50 |
!....................M8?????????????????????M.... | .........................
|
|
|
51 |
!.....................:M????????????????????M ... | DP .........................
|
|
|
52 |
!.......................MO???????????????????M... | .........................
|
|
|
53 |
!........................ M??????????????????M~.. | .............................
|
|
|
54 |
!..........................MM????????????????IM.. | .............................
|
|
|
55 |
!............................MI????????????????M. v .............................
|
|
|
56 |
!.............................:MI???????????????M. .............................
|
|
|
57 |
!...............................DMI??????????????M$ .............................
|
|
|
58 |
!................................ MM??????????????DM ............................
|
|
|
59 |
!.................................. MMI????????????IM ...........................
|
|
|
60 |
!.....................................MMI????????????M ..........................
|
|
|
61 |
!........................................MM??????????IM..........................
|
|
|
62 |
!......................................... $MMI???????M .........................
|
|
|
63 |
!........................................... MMMMMMM7. ..........................
|
|
|
64 |
!............................................. .............................
|
|
|
65 |
!............................................. PMAX .............................
|
|
|
66 |
!............................................. .............................
|
|
|
67 |
!................................................................................
|
|
|
68 |
|
|
|
69 |
USE mo_f2kcli ! command line interface
|
|
|
70 |
USE netcdf_tools
|
|
|
71 |
|
|
|
72 |
|
|
|
73 |
|
|
|
74 |
implicit none
|
|
|
75 |
|
|
|
76 |
! VERSION
|
|
|
77 |
character(len=*), parameter :: VERSION = '1.4'
|
|
|
78 |
! ... INTERNAL parameter
|
|
|
79 |
integer, parameter :: str_short = 30 ! length of short strings
|
|
|
80 |
integer, parameter :: str_long = 100 ! length of long strings
|
|
|
81 |
integer, parameter :: str_vlong = 500 ! length of very long strings
|
|
|
82 |
integer, parameter :: iou = 21 ! I/O unit
|
|
|
83 |
|
|
|
84 |
! FOR COMMAND LINE
|
|
|
85 |
CHARACTER(LEN=256) :: EXE ! program name
|
|
|
86 |
CHARACTER(LEN=80) :: CMD ! argument
|
|
|
87 |
INTEGER :: NARG ! number of arguments
|
|
|
88 |
|
|
|
89 |
! VARIABLES
|
|
|
90 |
INTEGER :: status ! status flag
|
|
|
91 |
LOGICAL :: l_invert = .FALSE. ! invert axis
|
|
|
92 |
LOGICAL :: l_Pa = .FALSE. ! convert units Pa->hPa
|
|
|
93 |
|
|
|
94 |
! Set tropopause thresholds
|
|
|
95 |
real :: tropo_pv,tropo_th
|
|
|
96 |
! diabatically produced PV (Q-threshold), Feb 2014:
|
|
|
97 |
real :: q_th
|
|
|
98 |
! masking number
|
|
|
99 |
real, parameter :: mdv =-999.0
|
|
|
100 |
|
|
|
101 |
! DATA/VARIABLE
|
|
|
102 |
! - input data
|
|
|
103 |
real, allocatable, dimension(:) :: x_data, y_data,z_data,t_data
|
|
|
104 |
real, allocatable, dimension(:) :: ak,bk
|
|
|
105 |
real, allocatable, dimension(:,:,:) :: aps
|
|
|
106 |
real, allocatable, dimension(:,:,:,:) :: press
|
|
|
107 |
real, allocatable, dimension(:,:,:,:) :: q_in,pv_in,pt_in
|
|
|
108 |
real, allocatable, dimension(:,:,:,:) :: q,pv,pt
|
|
|
109 |
character(len=str_long) :: x_units,y_units, z_units,t_units
|
|
|
110 |
character(len=str_long) :: aps_units
|
|
|
111 |
! - final results
|
|
|
112 |
real, allocatable, dimension(:,:,:) :: fold ! fold yes/no
|
|
|
113 |
real, allocatable, dimension(:,:,:) :: sfold ! shallow fold yes/no
|
|
|
114 |
real, allocatable, dimension(:,:,:) :: mfold ! medium fold yes/no
|
|
|
115 |
real, allocatable, dimension(:,:,:) :: dfold ! deep fold yes/no
|
|
|
116 |
real, allocatable, dimension(:,:,:) :: tp ! tropopause
|
|
|
117 |
real, allocatable, dimension(:,:,:) :: dp ! folding depth (Pa)
|
|
|
118 |
real, allocatable, dimension(:,:,:) :: pmin ! folding depth (Pa)
|
|
|
119 |
real, allocatable, dimension(:,:,:) :: pmax ! folding depth (Pa)
|
|
|
120 |
real, allocatable, dimension(:,:,:,:) :: label ! folding classification
|
|
|
121 |
real, allocatable, dimension(:,:,:,:) :: label_out
|
|
|
122 |
|
|
|
123 |
! Auxiliary variables
|
|
|
124 |
integer :: stat
|
|
|
125 |
integer :: i,j,k,t
|
|
|
126 |
real :: lev(100)
|
|
|
127 |
|
|
|
128 |
|
|
|
129 |
!namelist...
|
|
|
130 |
character(len=str_long) file_input
|
|
|
131 |
character(len=str_long) file_output
|
|
|
132 |
character(len=str_short) x_name,y_name,z_name,t_name
|
|
|
133 |
character(len=str_short) aps_name,hyam_name, hybm_name
|
|
|
134 |
character(len=str_short) q_name, pv_name, pt_name
|
|
|
135 |
|
|
|
136 |
! DATA DIMENSION
|
|
|
137 |
integer :: nx,ny,nz,nt
|
|
|
138 |
|
|
|
139 |
|
|
|
140 |
! -----------------------------------------------------------------
|
|
|
141 |
! Initialisation paramaters
|
|
|
142 |
! -----------------------------------------------------------------
|
|
|
143 |
! Humidity criterion for diabatically produced PV, BOJAN: this is the threshold Q=0.1g/kg
|
|
|
144 |
q_th=0.0001
|
|
|
145 |
tropo_pv = 2.0
|
|
|
146 |
tropo_th = 380.0
|
|
|
147 |
|
|
|
148 |
! -----------------------------------------------------------------
|
|
|
149 |
! Read command line
|
|
|
150 |
! -----------------------------------------------------------------
|
|
|
151 |
narg = command_argument_count() ! number of arguments
|
|
|
152 |
call get_command_argument(0,exe) ! program name
|
|
|
153 |
!
|
|
|
154 |
if (narg > 1) then
|
|
|
155 |
write(*,*) 'command-line error: too many arguments !'
|
|
|
156 |
call usage(trim(exe))
|
|
|
157 |
stop
|
|
|
158 |
end if
|
|
|
159 |
!
|
|
|
160 |
if (narg == 0) then
|
|
|
161 |
call usage(trim(exe))
|
|
|
162 |
stop
|
|
|
163 |
end if
|
|
|
164 |
!
|
|
|
165 |
call get_command_argument(1,cmd)
|
|
|
166 |
|
|
|
167 |
! -----------------------------------------------------------------
|
|
|
168 |
! Read namelist file
|
|
|
169 |
! -----------------------------------------------------------------
|
|
|
170 |
|
|
|
171 |
CALL read_nml(status, iou, TRIM(CMD))
|
|
|
172 |
IF (status /= 0) STOP
|
|
|
173 |
|
|
|
174 |
! -----------------------------------------------------------------
|
|
|
175 |
! Inquire netcdf file
|
|
|
176 |
! -----------------------------------------------------------------
|
|
|
177 |
CALL inquire_file(TRIM(file_input), &
|
|
|
178 |
x_name, y_name, z_name,t_name, &
|
|
|
179 |
nx, ny, nz, nt )
|
|
|
180 |
! IF (status > 0) STOP ! ERROR
|
|
|
181 |
|
|
|
182 |
|
|
|
183 |
WRITE(*,*) '==========================================================='
|
|
|
184 |
WRITE(*,*) " FILE DIMENSION: "
|
|
|
185 |
WRITE(*,*) " X :", nx
|
|
|
186 |
WRITE(*,*) " Y :", ny
|
|
|
187 |
WRITE(*,*) " Z :", nz
|
|
|
188 |
WRITE(*,*) " T :", nt
|
|
|
189 |
WRITE(*,*) '==========================================================='
|
|
|
190 |
|
|
|
191 |
! -----------------------------------------------------------------
|
|
|
192 |
! allocate data file
|
|
|
193 |
! -----------------------------------------------------------------
|
|
|
194 |
ALLOCATE(x_data(nx))
|
|
|
195 |
ALLOCATE(y_data(ny))
|
|
|
196 |
ALLOCATE(z_data(nz))
|
|
|
197 |
ALLOCATE(t_data(nt))
|
|
|
198 |
|
|
|
199 |
ALLOCATE(aps(nx,ny,nt))
|
|
|
200 |
ALLOCATE(ak(nz))
|
|
|
201 |
ALLOCATE(bk(nz))
|
|
|
202 |
ALLOCATE(press(nx,ny,nz,nt))
|
|
|
203 |
!
|
|
|
204 |
ALLOCATE(q(nx,ny,nz,nt))
|
|
|
205 |
ALLOCATE(pt(nx,ny,nz,nt))
|
|
|
206 |
ALLOCATE(pv(nx,ny,nz,nt))
|
|
|
207 |
ALLOCATE(q_in(nx,ny,nz,nt))
|
|
|
208 |
ALLOCATE(pt_in(nx,ny,nz,nt))
|
|
|
209 |
ALLOCATE(pv_in(nx,ny,nz,nt))
|
|
|
210 |
|
|
|
211 |
ALLOCATE(label(nx,ny,nz,nt))
|
|
|
212 |
ALLOCATE(fold(nx,ny,nt))
|
|
|
213 |
ALLOCATE(sfold(nx,ny,nt))
|
|
|
214 |
ALLOCATE(mfold(nx,ny,nt))
|
|
|
215 |
ALLOCATE(dfold(nx,ny,nt))
|
|
|
216 |
ALLOCATE(tp(nx,ny,nt))
|
|
|
217 |
ALLOCATE(dp(nx,ny,nt))
|
|
|
218 |
ALLOCATE(pmin(nx,ny,nt))
|
|
|
219 |
ALLOCATE(pmax(nx,ny,nt))
|
|
|
220 |
ALLOCATE(label_out(nx,ny,nz,nt))
|
|
|
221 |
|
|
|
222 |
! -----------------------------------------------------------------
|
|
|
223 |
! read netcdf file
|
|
|
224 |
! -----------------------------------------------------------------
|
|
|
225 |
call read_file(TRIM(file_input), TRIM(x_name), x_data)
|
|
|
226 |
call read_file(TRIM(file_input), TRIM(y_name), y_data)
|
|
|
227 |
call read_file(TRIM(file_input), TRIM(z_name), z_data)
|
|
|
228 |
call read_file(TRIM(file_input), TRIM(t_name), t_data)
|
|
|
229 |
|
|
|
230 |
call read_att(TRIM(file_input), TRIM(x_name), "units", x_units)
|
|
|
231 |
call read_att(TRIM(file_input), TRIM(y_name), "units", y_units)
|
|
|
232 |
call read_att(TRIM(file_input), TRIM(z_name), "units", z_units)
|
|
|
233 |
call read_att(TRIM(file_input), TRIM(t_name), "units", t_units)
|
|
|
234 |
|
|
|
235 |
! pressure data
|
|
|
236 |
call read_file(TRIM(file_input), TRIM(aps_name), aps)
|
|
|
237 |
call read_file(TRIM(file_input), TRIM(hyam_name), ak)
|
|
|
238 |
call read_file(TRIM(file_input), TRIM(hybm_name), bk)
|
|
|
239 |
! presure units
|
|
|
240 |
call read_att(TRIM(file_input), TRIM(aps_name), "units", aps_units)
|
|
|
241 |
|
|
|
242 |
! humidity
|
|
|
243 |
call read_file(TRIM(file_input), TRIM(q_name), q_in)
|
|
|
244 |
! potential temperature
|
|
|
245 |
call read_file(TRIM(file_input), TRIM(pt_name), pt_in)
|
|
|
246 |
! potential vorticity
|
|
|
247 |
call read_file(TRIM(file_input), TRIM(pv_name), pv_in)
|
|
|
248 |
|
|
|
249 |
! -----------------------------------------------------------------
|
|
|
250 |
! create pressure field and
|
|
|
251 |
! invert vertical axis (if needed)
|
|
|
252 |
! -----------------------------------------------------------------
|
|
|
253 |
! bottom => high pressure = 1
|
|
|
254 |
! top => low pressure = nz
|
|
|
255 |
if ((ak(1)+bk(1)*101325).lt.(ak(2)+bk(2)*101325)) then
|
|
|
256 |
l_invert=.TRUE.
|
|
|
257 |
write (*,*) 'Invert vertical axis'
|
|
|
258 |
endif
|
|
|
259 |
|
|
|
260 |
if (aps_units == "Pa") then
|
|
|
261 |
write (*,*) 'convert pressure from Pa to hPa'
|
|
|
262 |
l_Pa=.TRUE.
|
|
|
263 |
endif
|
|
|
264 |
|
|
|
265 |
if (l_invert) then
|
|
|
266 |
! to be inverted
|
|
|
267 |
do i=1,nx
|
|
|
268 |
do j=1,ny
|
|
|
269 |
do k=1,nz
|
|
|
270 |
do t=1,nt
|
|
|
271 |
if (l_Pa) then
|
|
|
272 |
press(i,j,nz-k+1,t)= (ak(k)+aps(i,j,t)*bk(k))/100.0 ! Pa -> hPa
|
|
|
273 |
else
|
|
|
274 |
press(i,j,nz-k+1,t)= (ak(k)+aps(i,j,t)*bk(k))
|
|
|
275 |
endif
|
|
|
276 |
q(i,j,nz-k+1,t) = q_in(i,j,k,t)
|
|
|
277 |
pt(i,j,nz-k+1,t)= pt_in(i,j,k,t)
|
|
|
278 |
pv(i,j,nz-k+1,t)= pv_in(i,j,k,t)
|
|
|
279 |
enddo
|
|
|
280 |
enddo
|
|
|
281 |
enddo
|
|
|
282 |
enddo
|
|
|
283 |
ELSE
|
|
|
284 |
! NOT to be inverted
|
|
|
285 |
do i=1,nx
|
|
|
286 |
do j=1,ny
|
|
|
287 |
do k=1,nz
|
|
|
288 |
do t=1,nt
|
|
|
289 |
if (l_Pa) then
|
|
|
290 |
press(i,j,k,t)= (ak(k)+aps(i,j,t)*bk(k))/100.0 ! Pa -> hPa
|
|
|
291 |
else
|
|
|
292 |
press(i,j,k,t)= ak(k)+aps(i,j,t)*bk(k)
|
|
|
293 |
endif
|
|
|
294 |
enddo
|
|
|
295 |
enddo
|
|
|
296 |
enddo
|
|
|
297 |
enddo
|
|
|
298 |
q(1:nx,1:ny,1:nz,1:nt) = q_in(1:nx,1:ny,1:nz,1:nt)
|
|
|
299 |
pt(1:nx,1:ny,1:nz,1:nt) = pt_in(1:nx,1:ny,1:nz,1:nt)
|
|
|
300 |
pv(1:nx,1:ny,1:nz,1:nt) = pv_in(1:nx,1:ny,1:nz,1:nt)
|
|
|
301 |
ENDIF
|
|
|
302 |
|
|
|
303 |
! -----------------------------------------------------------------
|
|
|
304 |
! perform calculation for each timestep
|
|
|
305 |
! -----------------------------------------------------------------
|
|
|
306 |
! TIME LOOP
|
|
|
307 |
do t=1,nt
|
|
|
308 |
|
|
|
309 |
|
|
|
310 |
write(*,*) "time step",t,"/", nt
|
|
|
311 |
! -------------------------------------------------------------------
|
|
|
312 |
! Part 1) Run clustering algorithm, get labels and interpolate to tropopause
|
|
|
313 |
! -------------------------------------------------------------------
|
|
|
314 |
|
|
|
315 |
! ------------------------------
|
|
|
316 |
! Running clustering algorithm
|
|
|
317 |
! ------------------------------
|
|
|
318 |
call tropopause (tp(:,:,t),label(:,:,:,t),press(:,:,:,t),pv(:,:,:,t),pt(:,:,:,t),q(:,:,:,t), &
|
|
|
319 |
mdv,y_data(:),nx,ny,nz,tropo_pv,tropo_th,q_th)
|
|
|
320 |
|
|
|
321 |
! -------------------------------------------------------------------
|
|
|
322 |
! Part 2) Identify tropopause folds from vertical cross-sections
|
|
|
323 |
! -----------------------------------------------------------------
|
|
|
324 |
|
|
|
325 |
! ------------------------------
|
|
|
326 |
! Identifying tropopause folds
|
|
|
327 |
! ------------------------------
|
|
|
328 |
call tropofold (label(:,:,:,t),press(:,:,:,t),pv(:,:,:,t),pt(:,:,:,t), &
|
|
|
329 |
fold(:,:,t), dp(:,:,t), pmin(:,:,t), pmax(:,:,t), sfold(:,:,t), mfold(:,:,t), dfold(:,:,t), &
|
|
|
330 |
mdv,nx,ny,nz,tropo_pv,tropo_th)
|
|
|
331 |
|
|
|
332 |
enddo
|
|
|
333 |
! END OF TIME LOOP
|
|
|
334 |
|
|
|
335 |
! -----------------------------------------------------------------
|
|
|
336 |
! re-invert vertical axis (if needed)
|
|
|
337 |
! -----------------------------------------------------------------
|
|
|
338 |
if (l_invert) then
|
|
|
339 |
! to be inverted
|
|
|
340 |
do i=1,nx
|
|
|
341 |
do j=1,ny
|
|
|
342 |
do k=1,nz
|
|
|
343 |
do t=1,nt
|
|
|
344 |
label_out(i,j,nz-k+1,t) = label(i,j,k,t)
|
|
|
345 |
enddo
|
|
|
346 |
enddo
|
|
|
347 |
enddo
|
|
|
348 |
enddo
|
|
|
349 |
ELSE
|
|
|
350 |
! NOT to be inverted
|
|
|
351 |
label_out(1:nx,1:ny,1:nz,1:nt) = label(1:nx,1:ny,1:nz,1:nt)
|
|
|
352 |
ENDIF
|
|
|
353 |
|
|
|
354 |
! -----------------------------------------------------------------
|
|
|
355 |
! write netcdf files
|
|
|
356 |
! -----------------------------------------------------------------
|
|
|
357 |
call nc_dump (TRIM(file_output), x_data, y_data, z_data, t_data, &
|
|
|
358 |
x_units,y_units, z_units,t_units, aps, ak, bk, label_out, fold, &
|
|
|
359 |
sfold, mfold, dfold, tp, dp, pmin, pmax)
|
|
|
360 |
|
|
|
361 |
|
|
|
362 |
! -----------------------------------------------------------------
|
|
|
363 |
! deallocate data file
|
|
|
364 |
! -----------------------------------------------------------------
|
|
|
365 |
DEALLOCATE(x_data)
|
|
|
366 |
DEALLOCATE(y_data)
|
|
|
367 |
DEALLOCATE(z_data)
|
|
|
368 |
DEALLOCATE(t_data)
|
|
|
369 |
|
|
|
370 |
DEALLOCATE(aps)
|
|
|
371 |
DEALLOCATE(ak)
|
|
|
372 |
DEALLOCATE(bk)
|
|
|
373 |
DEALLOCATE(press)
|
|
|
374 |
!
|
|
|
375 |
DEALLOCATE(q_in)
|
|
|
376 |
DEALLOCATE(pt_in)
|
|
|
377 |
DEALLOCATE(pv_in)
|
|
|
378 |
DEALLOCATE(q)
|
|
|
379 |
DEALLOCATE(pt)
|
|
|
380 |
DEALLOCATE(pv)
|
|
|
381 |
|
|
|
382 |
DEALLOCATE(label)
|
|
|
383 |
DEALLOCATE(label_out)
|
|
|
384 |
DEALLOCATE(fold)
|
|
|
385 |
DEALLOCATE(sfold)
|
|
|
386 |
DEALLOCATE(mfold)
|
|
|
387 |
DEALLOCATE(dfold)
|
|
|
388 |
DEALLOCATE(tp)
|
|
|
389 |
DEALLOCATE(dp)
|
|
|
390 |
DEALLOCATE(pmin)
|
|
|
391 |
DEALLOCATE(pmax)
|
|
|
392 |
|
|
|
393 |
CONTAINS
|
|
|
394 |
|
|
|
395 |
! ------------------------------------------------------------------------
|
|
|
396 |
SUBROUTINE USAGE(EXE)
|
|
|
397 |
CHARACTER (LEN=*) :: EXE
|
|
|
398 |
WRITE(*,*) '--------------------------------------------'
|
|
|
399 |
WRITE(*,*) '3d_labelling_and_fold_id, version: ',VERSION
|
|
|
400 |
WRITE(*,*) 'Authors: Andrea Pozzer, MPICH '
|
|
|
401 |
WRITE(*,*) ' Bojan Skerlak, ETH '
|
|
|
402 |
WRITE(*,*) ' Michael Sprenger, ETH '
|
|
|
403 |
WRITE(*,*) '--------------------------------------------'
|
|
|
404 |
WRITE(*,*) 'Usage: '//TRIM(EXE)//' <namelist-file>'
|
|
|
405 |
WRITE(*,*) 'See README.txt or EMAC.nml for example'
|
|
|
406 |
WRITE(*,*) '--------------------------------------------'
|
|
|
407 |
END SUBROUTINE USAGE
|
|
|
408 |
! ------------------------------------------------------------------------
|
|
|
409 |
|
|
|
410 |
! ------------------------------------------------------------------------
|
|
|
411 |
SUBROUTINE read_nml(status, iou, fname)
|
|
|
412 |
|
|
|
413 |
IMPLICIT NONE
|
|
|
414 |
|
|
|
415 |
! I/O
|
|
|
416 |
INTEGER, INTENT(OUT) :: status
|
|
|
417 |
INTEGER, INTENT(IN) :: iou
|
|
|
418 |
CHARACTER(LEN=*), INTENT(IN) :: fname
|
|
|
419 |
|
|
|
420 |
! LOCAL
|
|
|
421 |
CHARACTER(LEN=*), PARAMETER :: substr = 'read_nml'
|
|
|
422 |
LOGICAL :: lex ! file exists
|
|
|
423 |
INTEGER :: fstat ! file status
|
|
|
424 |
|
|
|
425 |
NAMELIST /CTRL/ file_input,file_output,X_name,Y_name,Z_name,T_name, &
|
|
|
426 |
APS_name,HYAM_name,HYBM_name,Q_name,PV_name,PT_name
|
|
|
427 |
|
|
|
428 |
status = 1 ! ERROR
|
|
|
429 |
|
|
|
430 |
WRITE(*,*) '==========================================================='
|
|
|
431 |
|
|
|
432 |
! CHECK IF FILE EXISTS
|
|
|
433 |
INQUIRE(file=TRIM(fname), exist=lex)
|
|
|
434 |
IF (.NOT.lex) THEN
|
|
|
435 |
WRITE(*,*) substr,': FILE DOES NOT EXIST (',TRIM(fname),')'
|
|
|
436 |
status = 1
|
|
|
437 |
RETURN
|
|
|
438 |
END IF
|
|
|
439 |
|
|
|
440 |
! OPEN FILE
|
|
|
441 |
OPEN(iou,file=TRIM(fname))
|
|
|
442 |
! READ NEMELIST
|
|
|
443 |
WRITE(*,*) 'READING NAMELIST ''CTRL'''//&
|
|
|
444 |
&' FROM '''//TRIM(fname),''' (unit ',iou,') ...'
|
|
|
445 |
!
|
|
|
446 |
READ(iou, NML=CTRL, IOSTAT=fstat)
|
|
|
447 |
!
|
|
|
448 |
IF (fstat /= 0) THEN
|
|
|
449 |
WRITE(*,*) substr,': READ ERROR IN NAMELIST ''CTRL'' (',TRIM(fname),')'
|
|
|
450 |
status = 3 ! READ ERROR IN NAMELIST
|
|
|
451 |
RETURN
|
|
|
452 |
END IF
|
|
|
453 |
|
|
|
454 |
WRITE(*,*) ' FILE INPUT : ', TRIM(file_input)
|
|
|
455 |
WRITE(*,*) ' FILE OUTPUT : ', TRIM(file_output)
|
|
|
456 |
WRITE(*,*) ' X_NAME : ', TRIM(x_name)
|
|
|
457 |
WRITE(*,*) ' Y_NAME : ', TRIM(y_name)
|
|
|
458 |
WRITE(*,*) ' Z_NAME : ', TRIM(z_name)
|
|
|
459 |
WRITE(*,*) ' T_NAME : ', TRIM(t_name)
|
|
|
460 |
WRITE(*,*) ' APS_NAME : ', TRIM(aps_name)
|
|
|
461 |
WRITE(*,*) ' HYAM_NAME : ', TRIM(hyam_name)
|
|
|
462 |
WRITE(*,*) ' HYBM_NAME : ', TRIM(hybm_name)
|
|
|
463 |
WRITE(*,*) ' Q_NAME : ', TRIM(q_name)
|
|
|
464 |
WRITE(*,*) ' PV_NAME : ', TRIM(pv_name)
|
|
|
465 |
WRITE(*,*) ' PT_NAME : ', TRIM(pt_name)
|
|
|
466 |
|
|
|
467 |
! CLOSE FILE
|
|
|
468 |
CLOSE(iou)
|
|
|
469 |
|
|
|
470 |
WRITE(*,*) '==========================================================='
|
|
|
471 |
|
|
|
472 |
status = 0
|
|
|
473 |
|
|
|
474 |
END SUBROUTINE read_nml
|
|
|
475 |
! ------------------------------------------------------------------------
|
|
|
476 |
|
|
|
477 |
|
|
|
478 |
end program clustering
|
|
|
479 |
|
|
|
480 |
|
|
|
481 |
! *******************************************************************************
|
|
|
482 |
! * Subroutine Section *
|
|
|
483 |
! *******************************************************************************
|
|
|
484 |
|
|
|
485 |
|
|
|
486 |
! --------------------------------------------------------
|
|
|
487 |
! Get grid points which are connected to a grid point: BOJAN: this is the 3D-connectedness criterion of the 3D-labelling algorithm
|
|
|
488 |
! --------------------------------------------------------
|
|
|
489 |
|
|
|
490 |
SUBROUTINE connect (outar,label,inx,iny,inz,dirh,dird,diru,nx,ny,nz)
|
|
|
491 |
|
|
|
492 |
! The input array <outar(nx,ny,nz)> is 0/1 labeled. Given a starting point
|
|
|
493 |
! <inx,iny,inz> where the label is 1, find all points with label 1 which are
|
|
|
494 |
! connected and attribute to all these points the label <label>.
|
|
|
495 |
! The 3D arrays <dir*(nx,ny,nz)> specifie whether an expansion is allowed
|
|
|
496 |
! dirh=1 where horizontal propagation is allowed,
|
|
|
497 |
! dird=1 where vertical propagation is allowed downward,
|
|
|
498 |
! diru=1 where vertical propagation is allowed upward, and
|
|
|
499 |
! a value of dir*=0 prohibits the respective propagation
|
|
|
500 |
|
|
|
501 |
! Declaration of subroutine parameters
|
|
|
502 |
integer :: nx,ny,nz
|
|
|
503 |
integer :: inx,iny,inz
|
|
|
504 |
integer :: outar(nx,ny,nz)
|
|
|
505 |
integer :: label
|
|
|
506 |
integer :: dirh(nx,ny,nz)
|
|
|
507 |
integer :: diru(nx,ny,nz)
|
|
|
508 |
integer :: dird(nx,ny,nz)
|
|
|
509 |
|
|
|
510 |
! Auxiliary variables
|
|
|
511 |
integer :: il,ir,jb,jf,ku,kd,im,jm,km
|
|
|
512 |
integer :: i,j,k
|
|
|
513 |
integer :: stack
|
|
|
514 |
integer :: indx(nx*ny*nz),indy(nx*ny*nz),indz(nx*ny*nz)
|
|
|
515 |
|
|
|
516 |
! Push the starting elements on the stack
|
|
|
517 |
stack=1
|
|
|
518 |
indx(stack)=inx
|
|
|
519 |
indy(stack)=iny
|
|
|
520 |
indz(stack)=inz
|
|
|
521 |
outar(inx,iny,inz)=label
|
|
|
522 |
|
|
|
523 |
! Define the indices of the neighbouring elements
|
|
|
524 |
100 continue
|
|
|
525 |
|
|
|
526 |
il=indx(stack)-1
|
|
|
527 |
if (il.lt.1) il=nx
|
|
|
528 |
ir=indx(stack)+1
|
|
|
529 |
if (ir.gt.nx) ir=1
|
|
|
530 |
jb=indy(stack)-1
|
|
|
531 |
if (jb.lt.1) jb=1
|
|
|
532 |
jf=indy(stack)+1
|
|
|
533 |
if (jf.gt.ny) jf=ny
|
|
|
534 |
ku=indz(stack)+1
|
|
|
535 |
if (ku.gt.nz) ku=nz
|
|
|
536 |
kd=indz(stack)-1
|
|
|
537 |
if (kd.lt.1) kd=1
|
|
|
538 |
im=indx(stack)
|
|
|
539 |
jm=indy(stack)
|
|
|
540 |
km=indz(stack)
|
|
|
541 |
stack=stack-1
|
|
|
542 |
|
|
|
543 |
! Check for index overflow
|
|
|
544 |
if (stack.ge.(nx*ny*nz-nx)) then
|
|
|
545 |
print*,'Stack overflow while clustering...'
|
|
|
546 |
stop
|
|
|
547 |
endif
|
|
|
548 |
|
|
|
549 |
! Mark all connected elements (build up the stack)
|
|
|
550 |
! Mark level below/above if dird/diru and dirh allowed
|
|
|
551 |
! BOJAN: U=up, D=down, R=right, L=left, M=mid, F=front, B=back (positions in a 27-element cube of nearest neighbours)
|
|
|
552 |
|
|
|
553 |
if ((dirh(im,jm,km).eq.1).and.(dird(im,jm,km).eq.1)) then
|
|
|
554 |
! below:
|
|
|
555 |
if ( outar(il,jf,kd).eq.1) then
|
|
|
556 |
outar(il,jf,kd)=label
|
|
|
557 |
stack=stack+1
|
|
|
558 |
indx(stack)=il
|
|
|
559 |
indy(stack)=jf
|
|
|
560 |
indz(stack)=kd
|
|
|
561 |
endif
|
|
|
562 |
if ( outar(im,jf,kd).eq.1) then
|
|
|
563 |
outar(im,jf,kd)=label
|
|
|
564 |
stack=stack+1
|
|
|
565 |
indx(stack)=im
|
|
|
566 |
indy(stack)=jf
|
|
|
567 |
indz(stack)=kd
|
|
|
568 |
endif
|
|
|
569 |
if ( outar(ir,jf,kd).eq.1) then
|
|
|
570 |
outar(ir,jf,kd)=label
|
|
|
571 |
stack=stack+1
|
|
|
572 |
indx(stack)=ir
|
|
|
573 |
indy(stack)=jf
|
|
|
574 |
indz(stack)=kd
|
|
|
575 |
endif
|
|
|
576 |
if (outar(il,jm,kd).eq.1) then
|
|
|
577 |
outar(il,jm,kd)=label
|
|
|
578 |
stack=stack+1
|
|
|
579 |
indx(stack)=il
|
|
|
580 |
indy(stack)=jm
|
|
|
581 |
indz(stack)=kd
|
|
|
582 |
endif
|
|
|
583 |
if (outar(ir,jm,kd).eq.1) then
|
|
|
584 |
outar(ir,jm,kd)=label
|
|
|
585 |
stack=stack+1
|
|
|
586 |
indx(stack)=ir
|
|
|
587 |
indy(stack)=jm
|
|
|
588 |
indz(stack)=kd
|
|
|
589 |
endif
|
|
|
590 |
if (outar(il,jb,kd).eq.1) then
|
|
|
591 |
outar(il,jb,kd)=label
|
|
|
592 |
stack=stack+1
|
|
|
593 |
indx(stack)=il
|
|
|
594 |
indy(stack)=jb
|
|
|
595 |
indz(stack)=kd
|
|
|
596 |
endif
|
|
|
597 |
if (outar(im,jb,kd).eq.1) then
|
|
|
598 |
outar(im,jb,kd)=label
|
|
|
599 |
stack=stack+1
|
|
|
600 |
indx(stack)=im
|
|
|
601 |
indy(stack)=jb
|
|
|
602 |
indz(stack)=kd
|
|
|
603 |
endif
|
|
|
604 |
if (outar(ir,jb,kd).eq.1) then
|
|
|
605 |
outar(ir,jb,kd)=label
|
|
|
606 |
stack=stack+1
|
|
|
607 |
indx(stack)=ir
|
|
|
608 |
indy(stack)=jb
|
|
|
609 |
indz(stack)=kd
|
|
|
610 |
endif
|
|
|
611 |
endif
|
|
|
612 |
if ((dirh(im,jm,km).eq.1).and.(diru(im,jm,km).eq.1)) then
|
|
|
613 |
if (outar(il,jf,ku).eq.1) then
|
|
|
614 |
outar(il,jf,ku)=label
|
|
|
615 |
stack=stack+1
|
|
|
616 |
indx(stack)=il
|
|
|
617 |
indy(stack)=jf
|
|
|
618 |
indz(stack)=ku
|
|
|
619 |
endif
|
|
|
620 |
if (outar(im,jf,ku).eq.1) then
|
|
|
621 |
outar(im,jf,ku)=label
|
|
|
622 |
stack=stack+1
|
|
|
623 |
indx(stack)=im
|
|
|
624 |
indy(stack)=jf
|
|
|
625 |
indz(stack)=ku
|
|
|
626 |
endif
|
|
|
627 |
if (outar(ir,jf,ku).eq.1) then
|
|
|
628 |
outar(ir,jf,ku)=label
|
|
|
629 |
stack=stack+1
|
|
|
630 |
indx(stack)=ir
|
|
|
631 |
indy(stack)=jf
|
|
|
632 |
indz(stack)=ku
|
|
|
633 |
endif
|
|
|
634 |
if (outar(il,jm,ku).eq.1) then
|
|
|
635 |
outar(il,jm,ku)=label
|
|
|
636 |
stack=stack+1
|
|
|
637 |
indx(stack)=il
|
|
|
638 |
indy(stack)=jm
|
|
|
639 |
indz(stack)=ku
|
|
|
640 |
endif
|
|
|
641 |
if (outar(ir,jm,ku).eq.1) then
|
|
|
642 |
outar(ir,jm,ku)=label
|
|
|
643 |
stack=stack+1
|
|
|
644 |
indx(stack)=ir
|
|
|
645 |
indy(stack)=jm
|
|
|
646 |
indz(stack)=ku
|
|
|
647 |
endif
|
|
|
648 |
if (outar(il,jb,ku).eq.1) then
|
|
|
649 |
outar(il,jb,ku)=label
|
|
|
650 |
stack=stack+1
|
|
|
651 |
indx(stack)=il
|
|
|
652 |
indy(stack)=jb
|
|
|
653 |
indz(stack)=ku
|
|
|
654 |
endif
|
|
|
655 |
if (outar(im,jb,ku).eq.1) then
|
|
|
656 |
outar(im,jb,ku)=label
|
|
|
657 |
stack=stack+1
|
|
|
658 |
indx(stack)=im
|
|
|
659 |
indy(stack)=jb
|
|
|
660 |
indz(stack)=ku
|
|
|
661 |
endif
|
|
|
662 |
if (outar(ir,jb,ku).eq.1) then
|
|
|
663 |
outar(ir,jb,ku)=label
|
|
|
664 |
stack=stack+1
|
|
|
665 |
indx(stack)=ir
|
|
|
666 |
indy(stack)=jb
|
|
|
667 |
indz(stack)=ku
|
|
|
668 |
endif
|
|
|
669 |
endif
|
|
|
670 |
! Mark level directly below and above if allowed (dird/u=1)
|
|
|
671 |
if (dird(im,jm,km).eq.1) then
|
|
|
672 |
if (outar(im,jm,kd).eq.1) then
|
|
|
673 |
outar(im,jm,kd)=label
|
|
|
674 |
stack=stack+1
|
|
|
675 |
indx(stack)=im
|
|
|
676 |
indy(stack)=jm
|
|
|
677 |
indz(stack)=kd
|
|
|
678 |
endif
|
|
|
679 |
endif
|
|
|
680 |
if (diru(im,jm,km).eq.1) then
|
|
|
681 |
if (outar(im,jm,ku).eq.1) then
|
|
|
682 |
outar(im,jm,ku)=label
|
|
|
683 |
stack=stack+1
|
|
|
684 |
indx(stack)=im
|
|
|
685 |
indy(stack)=jm
|
|
|
686 |
indz(stack)=ku
|
|
|
687 |
endif
|
|
|
688 |
endif
|
|
|
689 |
! Mark other points on same level if allowed (dirh=1)
|
|
|
690 |
if (dirh(im,jm,km).eq.1) then
|
|
|
691 |
if (outar(il,jf,km).eq.1) then
|
|
|
692 |
outar(il,jf,km)=label
|
|
|
693 |
stack=stack+1
|
|
|
694 |
indx(stack)=il
|
|
|
695 |
indy(stack)=jf
|
|
|
696 |
indz(stack)=km
|
|
|
697 |
endif
|
|
|
698 |
if (outar(im,jf,km).eq.1) then
|
|
|
699 |
outar(im,jf,km)=label
|
|
|
700 |
stack=stack+1
|
|
|
701 |
indx(stack)=im
|
|
|
702 |
indy(stack)=jf
|
|
|
703 |
indz(stack)=km
|
|
|
704 |
endif
|
|
|
705 |
if (outar(ir,jf,km).eq.1) then
|
|
|
706 |
outar(ir,jf,km)=label
|
|
|
707 |
stack=stack+1
|
|
|
708 |
indx(stack)=ir
|
|
|
709 |
indy(stack)=jf
|
|
|
710 |
indz(stack)=km
|
|
|
711 |
endif
|
|
|
712 |
if (outar(il,jm,km).eq.1) then
|
|
|
713 |
outar(il,jm,km)=label
|
|
|
714 |
stack=stack+1
|
|
|
715 |
indx(stack)=il
|
|
|
716 |
indy(stack)=jm
|
|
|
717 |
indz(stack)=km
|
|
|
718 |
endif
|
|
|
719 |
if (outar(ir,jm,km).eq.1) then
|
|
|
720 |
outar(ir,jm,km)=label
|
|
|
721 |
stack=stack+1
|
|
|
722 |
indx(stack)=ir
|
|
|
723 |
indy(stack)=jm
|
|
|
724 |
indz(stack)=km
|
|
|
725 |
endif
|
|
|
726 |
if (outar(il,jb,km).eq.1) then
|
|
|
727 |
outar(il,jb,km)=label
|
|
|
728 |
stack=stack+1
|
|
|
729 |
indx(stack)=il
|
|
|
730 |
indy(stack)=jb
|
|
|
731 |
indz(stack)=km
|
|
|
732 |
endif
|
|
|
733 |
if (outar(im,jb,km).eq.1) then
|
|
|
734 |
outar(im,jb,km)=label
|
|
|
735 |
stack=stack+1
|
|
|
736 |
indx(stack)=im
|
|
|
737 |
indy(stack)=jb
|
|
|
738 |
indz(stack)=km
|
|
|
739 |
endif
|
|
|
740 |
if (outar(ir,jb,km).eq.1) then
|
|
|
741 |
outar(ir,jb,km)=label
|
|
|
742 |
stack=stack+1
|
|
|
743 |
indx(stack)=ir
|
|
|
744 |
indy(stack)=jb
|
|
|
745 |
indz(stack)=km
|
|
|
746 |
endif
|
|
|
747 |
endif
|
|
|
748 |
|
|
|
749 |
if (stack.gt.0) goto 100
|
|
|
750 |
|
|
|
751 |
! Exit point
|
|
|
752 |
700 continue
|
|
|
753 |
|
|
|
754 |
return
|
|
|
755 |
|
|
|
756 |
end subroutine connect
|
|
|
757 |
|
|
|
758 |
! ---------------------------------------------------------
|
|
|
759 |
! Calculate the height of the tropopause: BOJAN: although this subroutine is called tropopause for historical reasons (long story), it actually is the core of the 3D-labelling algorithm :-)
|
|
|
760 |
! ---------------------------------------------------------
|
|
|
761 |
|
|
|
762 |
SUBROUTINE tropopause (f2,f3,p3,pv3,th3,q3,mdv,xlat,nx,ny,nz,tropo_pv,tropo_th,q_th)
|
|
|
763 |
|
|
|
764 |
! Get the pressure height of the tropopause (f2) and
|
|
|
765 |
! cluster the atmosphere into labels 1-5 (f3).
|
|
|
766 |
! Given: 3d field: Pressure <p3(nx,ny,nz)>,
|
|
|
767 |
! potential vorticity <pv3(nx,ny,nz)>, and potential
|
|
|
768 |
! temperature <th3(nx,ny,nz)>. The parameters <nx,ny,nz>
|
|
|
769 |
! characterize the grid. The missing data value is <mdv>.
|
|
|
770 |
! Bojan July 2013: tropo_pv and tropo_th have to be specified!
|
|
|
771 |
|
|
|
772 |
implicit none
|
|
|
773 |
|
|
|
774 |
! Declaration of subroutine parameters
|
|
|
775 |
integer,intent(in) :: nx,ny,nz
|
|
|
776 |
real,intent(out) :: f3(nx,ny,nz)
|
|
|
777 |
real,intent(out) :: f2(nx,ny)
|
|
|
778 |
real,intent(in) :: p3(nx,ny,nz)
|
|
|
779 |
real,intent(in) :: q3(nx,ny,nz)
|
|
|
780 |
real,intent(inout) :: pv3(nx,ny,nz)
|
|
|
781 |
real,intent(in) :: th3(nx,ny,nz)
|
|
|
782 |
real,intent(in) :: xlat(ny)
|
|
|
783 |
integer :: dirsh(nx,ny,nz),dirsd(nx,ny,nz),dirsu(nx,ny,nz)
|
|
|
784 |
integer :: dirth(nx,ny,nz),dirtd(nx,ny,nz),dirtu(nx,ny,nz)
|
|
|
785 |
integer :: dirdh(nx,ny,nz),dirdd(nx,ny,nz),dirdu(nx,ny,nz)
|
|
|
786 |
integer :: dirfh(nx,ny,nz),dirfd(nx,ny,nz),dirfu(nx,ny,nz)
|
|
|
787 |
real :: mdv
|
|
|
788 |
|
|
|
789 |
! Set tropopause thresholds
|
|
|
790 |
real :: tropo_pv,tropo_th
|
|
|
791 |
! diabatically produced PV (Q-threshold), Feb 2014:
|
|
|
792 |
real :: q_th
|
|
|
793 |
|
|
|
794 |
! Set permission for expansion of stratospheric label. A horizontal
|
|
|
795 |
! expansion is forbidden if less than <forbid_h> of the grid column
|
|
|
796 |
! below a grid point is tropospheric
|
|
|
797 |
real :: forbid_h
|
|
|
798 |
parameter (forbid_h=0.5)
|
|
|
799 |
|
|
|
800 |
! Internal variables
|
|
|
801 |
integer :: i,j,k,l
|
|
|
802 |
real :: pv2(nx,ny)
|
|
|
803 |
real :: th2(nx,ny)
|
|
|
804 |
real :: pvn(nx,ny),pvs(nx,ny)
|
|
|
805 |
real :: lat
|
|
|
806 |
integer :: st(nx,ny,nz),tr(nx,ny,nz),de(nx,ny,nz),fi(nx,ny,nz)
|
|
|
807 |
real :: sign
|
|
|
808 |
real :: lev(nx,ny)
|
|
|
809 |
integer :: tschk,kabove,kbelow,vertpvchk
|
|
|
810 |
real :: below, total
|
|
|
811 |
|
|
|
812 |
! ----- STEP 1a ----- BOJAN: these steps should be the same as described in the appendix of my PhD thesis
|
|
|
813 |
|
|
|
814 |
! Mark all points above 50 hPa as stratospheric
|
|
|
815 |
do i=1,nx
|
|
|
816 |
do j=1,ny
|
|
|
817 |
do k=1,nz
|
|
|
818 |
if (xlat(j).ge.0.) then
|
|
|
819 |
sign=1.
|
|
|
820 |
else
|
|
|
821 |
sign=-1.
|
|
|
822 |
endif
|
|
|
823 |
if (p3(i,j,k).lt.50.) then
|
|
|
824 |
pv3(i,j,k)=sign*1.5*tropo_pv
|
|
|
825 |
endif
|
|
|
826 |
enddo
|
|
|
827 |
enddo
|
|
|
828 |
enddo
|
|
|
829 |
|
|
|
830 |
! Mark points with PV > PV threshold (e.g. 2 pvu) OR TH > TH threshold (e.g. 380 K)
|
|
|
831 |
do i=1,nx
|
|
|
832 |
do j=1,ny
|
|
|
833 |
do k=1,nz
|
|
|
834 |
if ((abs(pv3(i,j,k)).ge.tropo_pv).or.(th3(i,j,k).ge.tropo_th)) then
|
|
|
835 |
st(i,j,k)=1 ! 'stratosphere'
|
|
|
836 |
tr(i,j,k)=0 ! 'tropopsphere'
|
|
|
837 |
de(i,j,k)=1 ! additional array for de-treatment
|
|
|
838 |
else
|
|
|
839 |
st(i,j,k)=0
|
|
|
840 |
tr(i,j,k)=1
|
|
|
841 |
de(i,j,k)=0
|
|
|
842 |
endif
|
|
|
843 |
enddo
|
|
|
844 |
enddo
|
|
|
845 |
enddo
|
|
|
846 |
|
|
|
847 |
! ----- STEP 1b -----
|
|
|
848 |
|
|
|
849 |
! Set the expansion permissions for the individual labels (1=allowed, 0=forbidden)
|
|
|
850 |
! h=horizontal, v=vertical
|
|
|
851 |
do i=1,nx
|
|
|
852 |
do j=1,ny
|
|
|
853 |
do k=1,nz
|
|
|
854 |
|
|
|
855 |
! Set permissions for stratospheric (st) label: always vertical, never horizontal
|
|
|
856 |
dirsh(i,j,k) = 0
|
|
|
857 |
dirsd(i,j,k) = 1
|
|
|
858 |
dirsu(i,j,k) = 1
|
|
|
859 |
|
|
|
860 |
! Set permissions for fill up thing label (fi): always horizontal, only downward, never upward
|
|
|
861 |
dirfh(i,j,k) = 1
|
|
|
862 |
dirfd(i,j,k) = 1
|
|
|
863 |
dirfu(i,j,k) = 0
|
|
|
864 |
|
|
|
865 |
! Set permissions for tropospheric label (can always propagate)
|
|
|
866 |
dirth(i,j,k) = 1
|
|
|
867 |
dirtd(i,j,k) = 1
|
|
|
868 |
dirtu(i,j,k) = 1
|
|
|
869 |
|
|
|
870 |
! Set permissions for deep fold (de) label: always horizontal, vertical only if air is dry
|
|
|
871 |
dirdh(i,j,k) = 1
|
|
|
872 |
! Allow vertical propagation of de label only for dry air (Q < Q-threshold)
|
|
|
873 |
if ( q3(i,j,k).le.q_th ) then
|
|
|
874 |
dirdd(i,j,k) = 1
|
|
|
875 |
dirdu(i,j,k) = 1
|
|
|
876 |
else
|
|
|
877 |
dirdd(i,j,k) = 0
|
|
|
878 |
dirdu(i,j,k) = 0
|
|
|
879 |
endif
|
|
|
880 |
!endif
|
|
|
881 |
|
|
|
882 |
enddo
|
|
|
883 |
enddo
|
|
|
884 |
enddo
|
|
|
885 |
|
|
|
886 |
! ----- STEP 2a -----
|
|
|
887 |
|
|
|
888 |
do i=1,nx
|
|
|
889 |
do j=1,ny
|
|
|
890 |
! Determine stratosphere by connectivity criterion: for all points (x/y) at highest
|
|
|
891 |
! model level with st=1, find all connected points with st=1 => they get label st=2
|
|
|
892 |
! this label can only propagate vertically and it is important that this happens before
|
|
|
893 |
! label st=5 can propagate from the bottom and overwrite st=1 with st=5
|
|
|
894 |
if ((st(i,j,nz).eq.1).and.(p3(i,j,nz).lt.500.)) then
|
|
|
895 |
call connect (st,2,i,j,nz,dirsh,dirsd,dirsu,nx,ny,nz)
|
|
|
896 |
endif
|
|
|
897 |
enddo
|
|
|
898 |
enddo
|
|
|
899 |
|
|
|
900 |
do i=1,nx
|
|
|
901 |
do j=1,ny
|
|
|
902 |
|
|
|
903 |
! ----- STEP 2b -----
|
|
|
904 |
|
|
|
905 |
! for de array, do the same thing. The big difference is that the horizontal propagation is
|
|
|
906 |
! always allowed. In the vertical direction, it is only allowed if RH<50%, ensuring that
|
|
|
907 |
! very deep tropopause folds, which are not connected by st are reached by de.
|
|
|
908 |
|
|
|
909 |
! This is needed because for some cases, stratospheric air (e.g. at the boundary of a deep-reaching
|
|
|
910 |
! streamer) is not identified as such (the label st=2 is not allowed to propagate
|
|
|
911 |
! horizontally). The de=2 label can always propagate horizontally and is in the end
|
|
|
912 |
! used to distinguish 'real' cutoffs (with st=1 and de=1) from
|
|
|
913 |
! deep reaching stratospheric air (with st=1 and de=2)
|
|
|
914 |
if ((de(i,j,nz).eq.1).and.(p3(i,j,nz).lt.500.)) then
|
|
|
915 |
call connect (de,2,i,j,nz,dirdh,dirdd,dirdu,nx,ny,nz)
|
|
|
916 |
endif
|
|
|
917 |
|
|
|
918 |
! ----- STEP 2c -----
|
|
|
919 |
|
|
|
920 |
! Determine troposphere analogeously: for all points (x/y) at lowest model
|
|
|
921 |
! level with tr=1 => all connected points get label tr=2
|
|
|
922 |
if ((tr(i,j,1).eq.1).and.(p3(i,j,1).gt.300.)) then
|
|
|
923 |
call connect (tr,2,i,j,1,dirth,dirtd,dirtu,nx,ny,nz)
|
|
|
924 |
endif
|
|
|
925 |
|
|
|
926 |
! ----- STEP 2d -----
|
|
|
927 |
|
|
|
928 |
! Determine surface-bound PV blobs (|PV|>2): for all points (x/y) at lowest
|
|
|
929 |
! model level with st=1 => all connected points get label st=5
|
|
|
930 |
! This label can always propagate! (use dirth,dirtv)
|
|
|
931 |
if ((st(i,j,1).eq.1).and.(p3(i,j,1).gt.300.)) then
|
|
|
932 |
call connect (st,5,i,j,1,dirth,dirtd,dirtu,nx,ny,nz)
|
|
|
933 |
endif
|
|
|
934 |
|
|
|
935 |
enddo
|
|
|
936 |
enddo
|
|
|
937 |
|
|
|
938 |
! ----- STEP 2e -----
|
|
|
939 |
|
|
|
940 |
! Remove tropospheric blobs in surface-bound PV blobs (tr=1 surrounded by st=5)
|
|
|
941 |
do i=1,nx
|
|
|
942 |
do j=1,ny
|
|
|
943 |
do k=1,nz
|
|
|
944 |
|
|
|
945 |
if (tr(i,j,k).eq.1) then
|
|
|
946 |
tschk=0
|
|
|
947 |
do l=k,nz
|
|
|
948 |
if (tschk.eq.0) then
|
|
|
949 |
! Check if air above is stratosphere (st=2) or top of atmosphere (nz)
|
|
|
950 |
if ( (st(i,j,l).eq.2) .or. (l.eq.nz) ) then
|
|
|
951 |
tschk=2
|
|
|
952 |
! Check if air above is surface-bound PV blob (st=5)
|
|
|
953 |
elseif (st(i,j,l).eq.5) then
|
|
|
954 |
tschk=5
|
|
|
955 |
endif
|
|
|
956 |
endif
|
|
|
957 |
enddo
|
|
|
958 |
! Mark this blob as tr=5
|
|
|
959 |
if (tschk.eq.5) then
|
|
|
960 |
call connect (tr,5,i,j,k,dirth,dirtd,dirtu,nx,ny,nz)
|
|
|
961 |
endif
|
|
|
962 |
endif
|
|
|
963 |
|
|
|
964 |
enddo
|
|
|
965 |
enddo
|
|
|
966 |
enddo
|
|
|
967 |
|
|
|
968 |
|
|
|
969 |
! ----- STEP 3 -----
|
|
|
970 |
|
|
|
971 |
! --------------------------------------------------------------------------------------------------
|
|
|
972 |
! Set the stratospheric and troposhperic labels (commented out version is 1=strat/0=everything else)
|
|
|
973 |
! --------------------------------------------------------------------------------------------------
|
|
|
974 |
do i=1,nx
|
|
|
975 |
do j=1,ny
|
|
|
976 |
do k=1,nz
|
|
|
977 |
f3(i,j,k)=0.
|
|
|
978 |
|
|
|
979 |
if (th3(i,j,k).gt.tropo_th) then
|
|
|
980 |
f3(i,j,k) = 2. ! stratosphere TH (2)
|
|
|
981 |
! f3(i,j,k) = 1. ! stratosphere TH (1)
|
|
|
982 |
else if (st(i,j,k).eq.2) then
|
|
|
983 |
f3(i,j,k) = 2. ! stratosphere PV (2)
|
|
|
984 |
! f3(i,j,k) = 1. ! stratosphere PV (1)
|
|
|
985 |
else if (tr(i,j,k).eq.2) then
|
|
|
986 |
f3(i,j,k) = 1. ! troposphere (1)
|
|
|
987 |
! f3(i,j,k) = 0. ! troposphere (0)
|
|
|
988 |
else if ((st(i,j,k).eq.1).and.(de(i,j,k).eq.1)) then
|
|
|
989 |
f3(i,j,k) = 3. ! strat cutoff (3)
|
|
|
990 |
! f3(i,j,k) = 0. ! strat cutoff (0)
|
|
|
991 |
else if ((st(i,j,k).eq.1).and.(de(i,j,k).eq.2)) then
|
|
|
992 |
f3(i,j,k) = 2. ! deep reaching stratospheric air (2)
|
|
|
993 |
! f3(i,j,k) = 1. ! deep reaching stratospheric air (1)
|
|
|
994 |
else if (tr(i,j,k).eq.1) then
|
|
|
995 |
f3(i,j,k) = 4. ! trop cutoff (4)
|
|
|
996 |
! f3(i,j,k) = 0. ! trop cutoff (0)
|
|
|
997 |
else if (st(i,j,k).eq.5) then
|
|
|
998 |
f3(i,j,k) = 5. ! surface-bound PV blob (5)
|
|
|
999 |
! f3(i,j,k) = 0. ! surface-bound PV blob (0)
|
|
|
1000 |
else if (tr(i,j,k).eq.5) then
|
|
|
1001 |
f3(i,j,k) = 1. ! tropospheric air within surface-bound PV blob (1)
|
|
|
1002 |
! f3(i,j,k) = 0. ! tropospheric air within surface-bound PV blob (0)
|
|
|
1003 |
endif
|
|
|
1004 |
|
|
|
1005 |
enddo
|
|
|
1006 |
enddo
|
|
|
1007 |
enddo
|
|
|
1008 |
|
|
|
1009 |
! --------------------------------------------------------------------------------------------------
|
|
|
1010 |
! --------------------------------------------------------------------------------------------------
|
|
|
1011 |
! ----- STEP 6 -----
|
|
|
1012 |
do i=1,nx
|
|
|
1013 |
do j=1,ny
|
|
|
1014 |
do k=1,nz
|
|
|
1015 |
fi(i,j,k)=0
|
|
|
1016 |
if (( f3(i,j,k).eq.2 ) .or. ( f3(i,j,k).eq.3 ).or.( f3(i,j,k).eq.5 )) then
|
|
|
1017 |
fi(i,j,k)=1 ! mark point as possible horizontal fill-up point
|
|
|
1018 |
endif
|
|
|
1019 |
end do
|
|
|
1020 |
end do
|
|
|
1021 |
end do
|
|
|
1022 |
! Horizontally fill up 2/3 and 2/5 patches (only aestetically important ;-))
|
|
|
1023 |
do i=1,nx
|
|
|
1024 |
do j=1,ny
|
|
|
1025 |
do k=1,nz
|
|
|
1026 |
if ( f3(i,j,k).eq.3 ) then ! label 3
|
|
|
1027 |
call connect (fi,3,i,j,k,dirfh,dirfd,dirfu,nx,ny,nz)
|
|
|
1028 |
endif
|
|
|
1029 |
if ( f3(i,j,k).eq.5 ) then ! label 3
|
|
|
1030 |
call connect (fi,5,i,j,k,dirfh,dirfd,dirfu,nx,ny,nz)
|
|
|
1031 |
endif
|
|
|
1032 |
enddo
|
|
|
1033 |
enddo
|
|
|
1034 |
enddo
|
|
|
1035 |
do i=1,nx
|
|
|
1036 |
do j=1,ny
|
|
|
1037 |
do k=1,nz
|
|
|
1038 |
if ( fi(i,j,k).eq.3 ) then ! label 3 that was horizontally filled up
|
|
|
1039 |
f3(i,j,k)=3.
|
|
|
1040 |
do l=1,k-1 ! vertically fill up points below if label is 2
|
|
|
1041 |
if ( f3(i,j,l).eq.2 ) then
|
|
|
1042 |
f3(i,j,l)=3.
|
|
|
1043 |
end if
|
|
|
1044 |
end do
|
|
|
1045 |
endif
|
|
|
1046 |
if ( fi(i,j,k).eq.5 ) then ! label 5 that was horizontally filled up
|
|
|
1047 |
f3(i,j,k)=5.
|
|
|
1048 |
endif
|
|
|
1049 |
enddo
|
|
|
1050 |
enddo
|
|
|
1051 |
enddo
|
|
|
1052 |
|
|
|
1053 |
! CHECK that all points obtained a label!
|
|
|
1054 |
do i=1,nx
|
|
|
1055 |
do j=1,ny
|
|
|
1056 |
do k=1,nz
|
|
|
1057 |
if ( f3(i,j,k).eq.0 ) then ! no label assigned
|
|
|
1058 |
print*,'Error: no label assigned'
|
|
|
1059 |
print*,'i,j,k = ',i,j,k
|
|
|
1060 |
print*,'p3(i,j,k) = ',p3(i,j,k)
|
|
|
1061 |
print*,'st(i,j,k) = ',st(i,j,k)
|
|
|
1062 |
print*,'tr(i,j,k) = ',tr(i,j,k)
|
|
|
1063 |
print*,'de(i,j,k) = ',de(i,j,k)
|
|
|
1064 |
stop
|
|
|
1065 |
endif
|
|
|
1066 |
enddo
|
|
|
1067 |
enddo
|
|
|
1068 |
enddo
|
|
|
1069 |
|
|
|
1070 |
! --------------------------------------------------------------------------------------------------
|
|
|
1071 |
! DONE with assigning labels, now interpolate field to tropopause pressure
|
|
|
1072 |
! --------------------------------------------------------------------------------------------------
|
|
|
1073 |
! BOJAN: here, we calculate the height of the tropopause (folds are ignored, see PhD thesis)
|
|
|
1074 |
! Remove stratospheric blobs in troposphere and tropospheric blobs in stratosphere
|
|
|
1075 |
! Brute force: Change PV values above and below tropo threshold
|
|
|
1076 |
do i=1,nx
|
|
|
1077 |
do j=1,ny
|
|
|
1078 |
do k=1,nz
|
|
|
1079 |
if (xlat(j).ge.0.) then
|
|
|
1080 |
sign=1.
|
|
|
1081 |
else
|
|
|
1082 |
sign=-1.
|
|
|
1083 |
endif
|
|
|
1084 |
|
|
|
1085 |
! Set PV of tropo. cutoff in stratosphere
|
|
|
1086 |
if ( (st(i,j,k).eq.0).and.(tr(i,j,k).eq.1) ) then
|
|
|
1087 |
pv3(i,j,k)=sign*1.5*tropo_pv
|
|
|
1088 |
! Set PV of exotic cases where PV < -2 pvu (PV > +2 pvu) in NH (SH) => would otherwise be detected as TP
|
|
|
1089 |
else if ( ((sign .gt. 0.).and.(pv3(i,j,k).lt.-tropo_pv).or.(sign .lt. 0.).and.(pv3(i,j,k).gt.tropo_pv)).and.(th3(i,j,k).gt.tropo_pv) ) then
|
|
|
1090 |
pv3(i,j,k)=sign*1.5*tropo_pv
|
|
|
1091 |
! Set PV of strat. cutoff in trop. or surface bound PV (only below 380K)
|
|
|
1092 |
!else if ( (((st(i,j,k).eq.1).and.(de(i,j,k).eq.1)).or.(st(i,j,k).eq.5)).and.(th3(i,j,k).le.tropo_th) ) then
|
|
|
1093 |
else if ( (f3(i,j,k).eq.3).or.(f3(i,j,k).eq.5) ) then
|
|
|
1094 |
pv3(i,j,k)=sign*0.5*tropo_pv
|
|
|
1095 |
endif
|
|
|
1096 |
enddo
|
|
|
1097 |
enddo
|
|
|
1098 |
enddo
|
|
|
1099 |
|
|
|
1100 |
! Calculate the height of the dynamical tropopause (NH:2pvu, SH:-2pvu)
|
|
|
1101 |
do i=1,nx
|
|
|
1102 |
do j=1,ny
|
|
|
1103 |
lev(i,j)=tropo_pv
|
|
|
1104 |
enddo
|
|
|
1105 |
enddo
|
|
|
1106 |
call thipo_lev(p3,pv3,lev,pvn,nx,ny,nz,mdv,0) ! for NH
|
|
|
1107 |
do i=1,nx
|
|
|
1108 |
do j=1,ny
|
|
|
1109 |
do k=1,nz
|
|
|
1110 |
pv3(i,j,k)=-pv3(i,j,k) ! change sign of PV
|
|
|
1111 |
enddo
|
|
|
1112 |
enddo
|
|
|
1113 |
enddo
|
|
|
1114 |
call thipo_lev(p3,pv3,lev,pvs,nx,ny,nz,mdv,0) ! for SH
|
|
|
1115 |
do i=1,nx
|
|
|
1116 |
do j=1,ny
|
|
|
1117 |
do k=1,nz
|
|
|
1118 |
pv3(i,j,k)=-pv3(i,j,k) ! change sign of PV back
|
|
|
1119 |
enddo
|
|
|
1120 |
enddo
|
|
|
1121 |
enddo
|
|
|
1122 |
do j=1,ny
|
|
|
1123 |
do i=1,nx
|
|
|
1124 |
if (xlat(j).ge.0.) then
|
|
|
1125 |
pv2(i,j)=pvn(i,j) ! NH
|
|
|
1126 |
else
|
|
|
1127 |
pv2(i,j)=pvs(i,j) ! SH
|
|
|
1128 |
endif
|
|
|
1129 |
enddo
|
|
|
1130 |
enddo
|
|
|
1131 |
|
|
|
1132 |
! Special treatment for PV towers (in whole column, |PV|>2)
|
|
|
1133 |
do i=1,nx
|
|
|
1134 |
do j=1,ny
|
|
|
1135 |
if (pv2(i,j).eq.mdv) then
|
|
|
1136 |
vertpvchk=0
|
|
|
1137 |
do k=1,nz
|
|
|
1138 |
if (abs(pv3(i,j,k)).ge.tropo_pv) vertpvchk=vertpvchk+1
|
|
|
1139 |
enddo
|
|
|
1140 |
if(vertpvchk.eq.nz) then
|
|
|
1141 |
!print*,'PV tower! Tropopause reaches ground: ',p3(i,j,1)
|
|
|
1142 |
pv2(i,j)=p3(i,j,1)
|
|
|
1143 |
endif
|
|
|
1144 |
endif
|
|
|
1145 |
enddo
|
|
|
1146 |
enddo
|
|
|
1147 |
|
|
|
1148 |
! Calculate the height of the 380 K surface
|
|
|
1149 |
do i=1,nx
|
|
|
1150 |
do j=1,ny
|
|
|
1151 |
lev(i,j)=tropo_th
|
|
|
1152 |
enddo
|
|
|
1153 |
enddo
|
|
|
1154 |
call thipo_lev(p3,th3,lev,th2,nx,ny,nz,mdv,0)
|
|
|
1155 |
|
|
|
1156 |
! Set unreasonable values to mdv
|
|
|
1157 |
do i=1,nx
|
|
|
1158 |
do j=1,ny
|
|
|
1159 |
if ((th2(i,j).lt.40.).or.(th2(i,j).gt.800.)) then ! Out of 'regular' limits => check label
|
|
|
1160 |
if (th2(i,j).ne.mdv) then
|
|
|
1161 |
!print*,'Warning: p(TP) (TH) out of limits:',th2(i,j)
|
|
|
1162 |
!print*,'i,j = ',i,j
|
|
|
1163 |
!print*,'LABEL profile:',f3(i,j,:)
|
|
|
1164 |
do k=1,nz
|
|
|
1165 |
if (p3(i,j,k).ge.th2(i,j)) then
|
|
|
1166 |
kbelow=k
|
|
|
1167 |
endif
|
|
|
1168 |
if (p3(i,j,nz+1-k).le.th2(i,j)) then
|
|
|
1169 |
kabove=nz+1-k
|
|
|
1170 |
endif
|
|
|
1171 |
enddo
|
|
|
1172 |
if ((f3(i,j,kbelow).eq.1).and.(f3(i,j,kabove).eq.2)) then
|
|
|
1173 |
!print*,'Point is okay!'
|
|
|
1174 |
else
|
|
|
1175 |
!print*,'Point is not okay => set to mdv'
|
|
|
1176 |
th2(i,j)=mdv
|
|
|
1177 |
endif
|
|
|
1178 |
endif
|
|
|
1179 |
endif
|
|
|
1180 |
if ((pv2(i,j).lt.40.).or.(pv2(i,j).gt.800.)) then ! Out of 'regular' limits => check label
|
|
|
1181 |
if (pv2(i,j).ne.mdv) then
|
|
|
1182 |
!print*,'Warning: p(TP) (PV) out of limits:',pv2(i,j)
|
|
|
1183 |
!print*,'i,j = ',i,j
|
|
|
1184 |
!print*,'LABEL profile:',f3(i,j,:)
|
|
|
1185 |
do k=1,nz
|
|
|
1186 |
if (p3(i,j,k).ge.pv2(i,j)) then
|
|
|
1187 |
kbelow=k
|
|
|
1188 |
endif
|
|
|
1189 |
if (p3(i,j,nz+1-k).le.pv2(i,j)) then
|
|
|
1190 |
kabove=nz+1-k
|
|
|
1191 |
endif
|
|
|
1192 |
enddo
|
|
|
1193 |
if ((f3(i,j,kbelow).eq.1).and.(f3(i,j,kabove).eq.2)) then
|
|
|
1194 |
!print*,'Point is okay!'
|
|
|
1195 |
elseif ((kabove.eq.1).and.(kbelow.eq.1)) then
|
|
|
1196 |
!print*,'PV tower! Tropopause reaches ground: ',pv2(i,j)
|
|
|
1197 |
else
|
|
|
1198 |
!print*,'Point is not okay => set to mdv'
|
|
|
1199 |
pv2(i,j)=mdv
|
|
|
1200 |
endif
|
|
|
1201 |
endif
|
|
|
1202 |
|
|
|
1203 |
endif
|
|
|
1204 |
enddo
|
|
|
1205 |
enddo
|
|
|
1206 |
|
|
|
1207 |
! --------------------------------------------------------------------------------------------------
|
|
|
1208 |
! Set the resulting tropopause height (maximum pressure: 380 K, 2 PVU)
|
|
|
1209 |
! --------------------------------------------------------------------------------------------------
|
|
|
1210 |
do i=1,nx
|
|
|
1211 |
do j=1,ny
|
|
|
1212 |
if ((th2(i,j).ne.mdv).and.(pv2(i,j).eq.mdv)) then
|
|
|
1213 |
f2(i,j)=th2(i,j)
|
|
|
1214 |
else if ((th2(i,j).eq.mdv).and.(pv2(i,j).ne.mdv)) then
|
|
|
1215 |
f2(i,j)=pv2(i,j)
|
|
|
1216 |
else if (th2(i,j).gt.pv2(i,j)) then
|
|
|
1217 |
f2(i,j)=th2(i,j)
|
|
|
1218 |
else
|
|
|
1219 |
f2(i,j)=pv2(i,j)
|
|
|
1220 |
endif
|
|
|
1221 |
enddo
|
|
|
1222 |
enddo
|
|
|
1223 |
! --------------------------------------------------------------------------------------------------
|
|
|
1224 |
! DONE
|
|
|
1225 |
! --------------------------------------------------------------------------------------------------
|
|
|
1226 |
|
|
|
1227 |
end subroutine tropopause
|
|
|
1228 |
|
|
|
1229 |
SUBROUTINE tropofold (label3,p3,pv3,th3,fold,dp, pmin, pmax, sfold,mfold,dfold,mdv,nx,ny,nz,tropo_pv,tropo_th)
|
|
|
1230 |
|
|
|
1231 |
! -------------------------------------------------------------------------------
|
|
|
1232 |
! Part 2) Identify tropopause folds (multiple crossings of the tropopause in one column)
|
|
|
1233 |
! -------------------------------------------------------------------------------
|
|
|
1234 |
|
|
|
1235 |
! Given the label field, determine whether a tropopause fold is present.
|
|
|
1236 |
! The parameters <nx,ny,nx> characterize the grid. The missing data value is <mdv>.
|
|
|
1237 |
|
|
|
1238 |
implicit none
|
|
|
1239 |
|
|
|
1240 |
! Declaration of subroutine parameters
|
|
|
1241 |
integer :: nx,ny,nz
|
|
|
1242 |
real :: mdv
|
|
|
1243 |
real :: label3(nx,ny,nz)
|
|
|
1244 |
real :: p3(nx,ny,nz)
|
|
|
1245 |
real :: pv3(nx,ny,nz)
|
|
|
1246 |
real :: th3(nx,ny,nz)
|
|
|
1247 |
real :: pre0,pre1,pre2,ok1,ok2
|
|
|
1248 |
real :: tropo_pv,tropo_th
|
|
|
1249 |
real :: fold(nx,ny)
|
|
|
1250 |
real :: dp(nx,ny)
|
|
|
1251 |
real :: pmin(nx,ny)
|
|
|
1252 |
real :: pmax(nx,ny)
|
|
|
1253 |
real :: sfold(nx,ny)
|
|
|
1254 |
real :: mfold(nx,ny)
|
|
|
1255 |
real :: dfold(nx,ny)
|
|
|
1256 |
|
|
|
1257 |
! Aux integers
|
|
|
1258 |
integer :: i,j,k,l,k0,k1,k2,ncross
|
|
|
1259 |
|
|
|
1260 |
! Adjust fields (integer labels, absolute value of PV)
|
|
|
1261 |
do i=1,nx
|
|
|
1262 |
do j=1,ny
|
|
|
1263 |
do k=1,nz
|
|
|
1264 |
pv3(i,j,k) = abs(pv3(i,j,k))
|
|
|
1265 |
enddo
|
|
|
1266 |
enddo
|
|
|
1267 |
enddo
|
|
|
1268 |
|
|
|
1269 |
! Definition of labels
|
|
|
1270 |
! 1 : troposphere
|
|
|
1271 |
! 2 : stratosphere
|
|
|
1272 |
! 3 : stratosperic cutoff
|
|
|
1273 |
! 4 : tropospheric cutoff
|
|
|
1274 |
! 5 : surface-bound PV blobs
|
|
|
1275 |
|
|
|
1276 |
do i=1,nx
|
|
|
1277 |
do j=1,ny
|
|
|
1278 |
k0 = 0
|
|
|
1279 |
k1 = 0
|
|
|
1280 |
k2 = 0
|
|
|
1281 |
ncross = 0
|
|
|
1282 |
|
|
|
1283 |
do k=nz-1,1,-1
|
|
|
1284 |
|
|
|
1285 |
! strat-trop transition
|
|
|
1286 |
if ( (label3(i,j,k+1).eq.2).and.(label3(i,j,k).eq.1) ) then
|
|
|
1287 |
|
|
|
1288 |
ncross = ncross + 1
|
|
|
1289 |
if ( k2.eq.0 ) then
|
|
|
1290 |
k2 = k
|
|
|
1291 |
elseif ( k1.ne.0 ) then
|
|
|
1292 |
k0 = k
|
|
|
1293 |
endif
|
|
|
1294 |
endif
|
|
|
1295 |
! special case for folds that are affected by the q-criterion (q<0.0001 only)
|
|
|
1296 |
! example (top to bottom): 222211122231111 need to recognize the 2-3-1 transition
|
|
|
1297 |
! transitions like 222222223111 should not be recognized, thus only if k2.ne.0
|
|
|
1298 |
if ( (label3(i,j,k+1).eq.2).and.(label3(i,j,k).eq.3) ) then
|
|
|
1299 |
|
|
|
1300 |
ncross = ncross + 1
|
|
|
1301 |
if (( k2.ne.0 ) .and. ( k1.ne.0 )) then
|
|
|
1302 |
k0 = k
|
|
|
1303 |
endif
|
|
|
1304 |
endif
|
|
|
1305 |
|
|
|
1306 |
! trop-strat transition
|
|
|
1307 |
if ( (label3(i,j,k+1).eq.1).and.(label3(i,j,k).eq.2) ) then
|
|
|
1308 |
|
|
|
1309 |
ncross = ncross + 1
|
|
|
1310 |
if ( ( k2.ne.0 ).and.( k1.eq.0 ) ) then
|
|
|
1311 |
k1 = k
|
|
|
1312 |
endif
|
|
|
1313 |
endif
|
|
|
1314 |
enddo
|
|
|
1315 |
|
|
|
1316 |
! Skip further steps if we have a single tropopause
|
|
|
1317 |
if ( ncross.le.2 ) goto 100
|
|
|
1318 |
! Get the exact pressures at the crossings
|
|
|
1319 |
pre0 = 0.
|
|
|
1320 |
pre1 = 0.
|
|
|
1321 |
pre2 = 0.
|
|
|
1322 |
! lowest point (0) => pre0
|
|
|
1323 |
ok1 = ( pv3(i,j,k0+1)-tropo_pv ) * ( pv3(i,j,k0)-tropo_pv )
|
|
|
1324 |
ok2 = ( th3(i,j,k0+1)-tropo_th ) * ( th3(i,j,k0)-tropo_th )
|
|
|
1325 |
|
|
|
1326 |
if ( ok1.le.0. ) then
|
|
|
1327 |
pre0 = p3(i,j,k0) + ( p3(i,j,k0+1) - p3(i,j,k0) ) * &
|
|
|
1328 |
( tropo_pv - pv3(i,j,k0) )/( pv3(i,j,k0+1) - pv3(i,j,k0) )
|
|
|
1329 |
elseif ( ok2.le.0. ) then
|
|
|
1330 |
pre0 = p3(i,j,k0) + ( p3(i,j,k0+1) - p3(i,j,k0) ) * &
|
|
|
1331 |
( tropo_th - th3(i,j,k0) )/( th3(i,j,k0+1) - th3(i,j,k0) )
|
|
|
1332 |
endif
|
|
|
1333 |
! middle point (1) => pre1
|
|
|
1334 |
ok1 = ( pv3(i,j,k1+1)-tropo_pv ) * ( pv3(i,j,k1)-tropo_pv )
|
|
|
1335 |
ok2 = ( th3(i,j,k1+1)-tropo_th ) * ( th3(i,j,k1)-tropo_th )
|
|
|
1336 |
if ( ok1.le.0. ) then
|
|
|
1337 |
pre1 = p3(i,j,k1) + ( p3(i,j,k1+1) - p3(i,j,k1) ) * &
|
|
|
1338 |
( tropo_pv - pv3(i,j,k1) )/( pv3(i,j,k1+1) - pv3(i,j,k1) )
|
|
|
1339 |
elseif ( ok2.le.0. ) then
|
|
|
1340 |
pre1 = p3(i,j,k1) + ( p3(i,j,k1+1) - p3(i,j,k1) ) * &
|
|
|
1341 |
( tropo_th - th3(i,j,k1) )/( th3(i,j,k1+1) - th3(i,j,k1) )
|
|
|
1342 |
endif
|
|
|
1343 |
! upper point (2) => pre2
|
|
|
1344 |
ok1 = ( pv3(i,j,k2+1)-tropo_pv ) * ( pv3(i,j,k2)-tropo_pv )
|
|
|
1345 |
ok2 = ( th3(i,j,k2+1)-tropo_th ) * ( th3(i,j,k2)-tropo_th )
|
|
|
1346 |
if ( ok1.le.0. ) then
|
|
|
1347 |
pre2 = p3(i,j,k2) + ( p3(i,j,k2+1) - p3(i,j,k2) ) * &
|
|
|
1348 |
( tropo_pv - pv3(i,j,k2) )/( pv3(i,j,k2+1) - pv3(i,j,k2) )
|
|
|
1349 |
elseif ( ok2.le.0. ) then
|
|
|
1350 |
pre2 = p3(i,j,k2) + ( p3(i,j,k2+1) - p3(i,j,k2) ) * &
|
|
|
1351 |
( tropo_th - th3(i,j,k2) )/( th3(i,j,k2+1) - th3(i,j,k2) )
|
|
|
1352 |
endif
|
|
|
1353 |
|
|
|
1354 |
! Decide whether all pressure values are ok
|
|
|
1355 |
if ( ( pre0.lt.p3(i,j,k0+1) ).or. &
|
|
|
1356 |
( pre0.gt.p3(i,j,k0 ) ).or. &
|
|
|
1357 |
( pre1.lt.p3(i,j,k1+1) ).or. &
|
|
|
1358 |
( pre1.gt.p3(i,j,k1 ) ).or. &
|
|
|
1359 |
( pre2.lt.p3(i,j,k2+1) ).or. &
|
|
|
1360 |
( pre2.gt.p3(i,j,k2 ) ) ) goto 100
|
|
|
1361 |
! Everything is fine: remember the fold
|
|
|
1362 |
fold(i,j) = 1.
|
|
|
1363 |
dp(i,j) = pre1 - pre2
|
|
|
1364 |
pmin(i,j) = pre2
|
|
|
1365 |
pmax(i,j) = pre0
|
|
|
1366 |
|
|
|
1367 |
! Exit point for loop
|
|
|
1368 |
100 continue
|
|
|
1369 |
|
|
|
1370 |
enddo ! y
|
|
|
1371 |
enddo ! x
|
|
|
1372 |
|
|
|
1373 |
where((fold .gt. 0.0) .and. (dp .ge. 50.0) .and. (dp .lt. 200.0))
|
|
|
1374 |
sfold=1.0
|
|
|
1375 |
elsewhere((fold .gt. 0.0) .and. (dp .ge. 200.0) .and. (dp .lt. 350.0))
|
|
|
1376 |
mfold=1.0
|
|
|
1377 |
elsewhere((fold .gt. 0.0) .and. (dp .ge. 350.0))
|
|
|
1378 |
dfold=1.0
|
|
|
1379 |
end where
|
|
|
1380 |
|
|
|
1381 |
end subroutine tropofold
|
|
|
1382 |
|
|
|
1383 |
! -----------------------------------------------------------
|
|
|
1384 |
! Interpolation onto theta surface (top -> down): BOJAN: not sure you have all the interpolation subroutines used here, contact us or feel free to use any other interpolation routine that does the job :-)
|
|
|
1385 |
! -----------------------------------------------------------
|
|
|
1386 |
subroutine thipo_lev(var3d,th3d,lev,var,nx,ny,nz,mdv,mode)
|
|
|
1387 |
|
|
|
1388 |
! Interpolates the 3d variable var3d on the isentropic surface
|
|
|
1389 |
! defined by lev. The interpolated field is returned as var.
|
|
|
1390 |
! th3d denotes the 3d theta array.
|
|
|
1391 |
! mode determines the way of vertical interpolation:
|
|
|
1392 |
! mode=0 is for linear interpolation
|
|
|
1393 |
! mode=1 is for logarithmic interpolation
|
|
|
1394 |
|
|
|
1395 |
integer :: nx,ny,nz,mode
|
|
|
1396 |
real :: lev(nx,ny),mdv
|
|
|
1397 |
real :: var3d(nx,ny,nz),th3d(nx,ny,nz),var(nx,ny)
|
|
|
1398 |
|
|
|
1399 |
integer :: i,j,k
|
|
|
1400 |
real :: kind
|
|
|
1401 |
real :: int3dm
|
|
|
1402 |
|
|
|
1403 |
do i=1,nx
|
|
|
1404 |
do j=1,ny
|
|
|
1405 |
|
|
|
1406 |
if (lev(i,j).ne.mdv) then
|
|
|
1407 |
|
|
|
1408 |
kind=0.
|
|
|
1409 |
do k=nz-1,1,-1
|
|
|
1410 |
if ((th3d(i,j,k).le.lev(i,j)).and.(th3d(i,j,k+1).ge.lev(i,j))) then
|
|
|
1411 |
kind=float(k)+(th3d(i,j,k)-lev(i,j))/(th3d(i,j,k)-th3d(i,j,k+1))
|
|
|
1412 |
goto 100
|
|
|
1413 |
endif
|
|
|
1414 |
enddo
|
|
|
1415 |
100 continue
|
|
|
1416 |
|
|
|
1417 |
if (kind.eq.0) then
|
|
|
1418 |
var(i,j)=mdv
|
|
|
1419 |
else
|
|
|
1420 |
var(i,j)=int3dm(var3d,nx,ny,nz,float(i),float(j),kind,mdv)
|
|
|
1421 |
endif
|
|
|
1422 |
|
|
|
1423 |
else
|
|
|
1424 |
|
|
|
1425 |
var(i,j)=mdv
|
|
|
1426 |
|
|
|
1427 |
endif
|
|
|
1428 |
|
|
|
1429 |
enddo
|
|
|
1430 |
enddo
|
|
|
1431 |
|
|
|
1432 |
return
|
|
|
1433 |
end subroutine thipo_lev
|
|
|
1434 |
|
|
|
1435 |
! -----------------------------------------------------------
|
|
|
1436 |
! Interpolation onto theta surface (bottom -> up)
|
|
|
1437 |
! -----------------------------------------------------------
|
|
|
1438 |
|
|
|
1439 |
subroutine pipo_lev(var3d,p3d,lev,var,nx,ny,nz,mdv,mode)
|
|
|
1440 |
|
|
|
1441 |
! Interpolates the 3d variable var3d on the pressure surface
|
|
|
1442 |
! defined by lev. The interpolated field is returned as var.
|
|
|
1443 |
! p3d denotes the 3d pressure array.
|
|
|
1444 |
! mode determines the way of vertical interpolation:
|
|
|
1445 |
! mode=0 is for linear interpolation
|
|
|
1446 |
! mode=1 is for logarithmic interpolation
|
|
|
1447 |
|
|
|
1448 |
integer :: nx,ny,nz,mode
|
|
|
1449 |
real :: lev(nx,ny),mdv
|
|
|
1450 |
real :: var3d(nx,ny,nz),p3d(nx,ny,nz),var(nx,ny)
|
|
|
1451 |
|
|
|
1452 |
integer :: i,j,k
|
|
|
1453 |
real :: kind
|
|
|
1454 |
real :: int3dm
|
|
|
1455 |
|
|
|
1456 |
do i=1,nx
|
|
|
1457 |
do j=1,ny
|
|
|
1458 |
|
|
|
1459 |
if (lev(i,j).ne.mdv) then
|
|
|
1460 |
|
|
|
1461 |
kind=0.
|
|
|
1462 |
do k=1,nz-1
|
|
|
1463 |
if ((p3d(i,j,k).ge.lev(i,j)).and.(p3d(i,j,k+1).le.lev(i,j))) then
|
|
|
1464 |
kind=float(k)+(p3d(i,j,k)-lev(i,j))/(p3d(i,j,k)-p3d(i,j,k+1))
|
|
|
1465 |
goto 100
|
|
|
1466 |
endif
|
|
|
1467 |
enddo
|
|
|
1468 |
100 continue
|
|
|
1469 |
|
|
|
1470 |
if (kind.eq.0.) then
|
|
|
1471 |
var(i,j)=mdv
|
|
|
1472 |
else
|
|
|
1473 |
var(i,j)=int3dm(var3d,nx,ny,nz,float(i),float(j),kind,mdv)
|
|
|
1474 |
endif
|
|
|
1475 |
|
|
|
1476 |
else
|
|
|
1477 |
|
|
|
1478 |
var(i,j)=mdv
|
|
|
1479 |
|
|
|
1480 |
endif
|
|
|
1481 |
|
|
|
1482 |
enddo
|
|
|
1483 |
enddo
|
|
|
1484 |
|
|
|
1485 |
return
|
|
|
1486 |
end subroutine pipo_lev
|
|
|
1487 |
|
|
|
1488 |
|
|
|
1489 |
! ------------------------------------------------------------------
|
|
|
1490 |
! Auxiliary routines
|
|
|
1491 |
! ------------------------------------------------------------------
|
|
|
1492 |
|