1 | !> @file data_output_3d.f90 |
---|
2 | !------------------------------------------------------------------------------! |
---|
3 | ! This file is part of the PALM model system. |
---|
4 | ! |
---|
5 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
6 | ! terms of the GNU General Public License as published by the Free Software |
---|
7 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
8 | ! version. |
---|
9 | ! |
---|
10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
13 | ! |
---|
14 | ! You should have received a copy of the GNU General Public License along with |
---|
15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
16 | ! |
---|
17 | ! Copyright 1997-2017 Leibniz Universitaet Hannover |
---|
18 | !------------------------------------------------------------------------------! |
---|
19 | ! |
---|
20 | ! Current revisions: |
---|
21 | ! ------------------ |
---|
22 | ! |
---|
23 | ! |
---|
24 | ! Former revisions: |
---|
25 | ! ----------------- |
---|
26 | ! $Id: data_output_3d.f90 2696 2017-12-14 17:12:51Z maronga $ |
---|
27 | ! Implementation of turbulence_closure_mod (TG) |
---|
28 | ! Implementation of chemistry module (FK) |
---|
29 | ! Set fill values at topography grid points or e.g. non-natural-type surface |
---|
30 | ! in case of LSM output (MS) |
---|
31 | ! |
---|
32 | ! 2512 2017-10-04 08:26:59Z raasch |
---|
33 | ! upper bounds of 3d output changed from nx+1,ny+1 to nx,ny |
---|
34 | ! no output of ghost layer data |
---|
35 | ! |
---|
36 | ! 2292 2017-06-20 09:51:42Z schwenkel |
---|
37 | ! Implementation of new microphysic scheme: cloud_scheme = 'morrison' |
---|
38 | ! includes two more prognostic equations for cloud drop concentration (nc) |
---|
39 | ! and cloud water content (qc). |
---|
40 | ! |
---|
41 | ! 2233 2017-05-30 18:08:54Z suehring |
---|
42 | ! |
---|
43 | ! 2232 2017-05-30 17:47:52Z suehring |
---|
44 | ! Adjustments to new topography concept |
---|
45 | ! |
---|
46 | ! 2209 2017-04-19 09:34:46Z kanani |
---|
47 | ! Added plant canopy model output |
---|
48 | ! |
---|
49 | ! 2031 2016-10-21 15:11:58Z knoop |
---|
50 | ! renamed variable rho to rho_ocean and rho_av to rho_ocean_av |
---|
51 | ! |
---|
52 | ! 2011 2016-09-19 17:29:57Z kanani |
---|
53 | ! Flag urban_surface is now defined in module control_parameters, |
---|
54 | ! changed prefix for urban surface model output to "usm_", |
---|
55 | ! introduced control parameter varnamelength for LEN of trimvar. |
---|
56 | ! |
---|
57 | ! 2007 2016-08-24 15:47:17Z kanani |
---|
58 | ! Added support for new urban surface model (temporary modifications of |
---|
59 | ! SELECT CASE ( ) necessary, see variable trimvar) |
---|
60 | ! |
---|
61 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
62 | ! Forced header and separation lines into 80 columns |
---|
63 | ! |
---|
64 | ! 1980 2016-07-29 15:51:57Z suehring |
---|
65 | ! Bugfix, in order to steer user-defined output, setting flag found explicitly |
---|
66 | ! to .F. |
---|
67 | ! |
---|
68 | ! 1976 2016-07-27 13:28:04Z maronga |
---|
69 | ! Output of radiation quantities is now done directly in the respective module |
---|
70 | ! |
---|
71 | ! 1972 2016-07-26 07:52:02Z maronga |
---|
72 | ! Output of land surface quantities is now done directly in the respective module. |
---|
73 | ! Unnecessary directive __parallel removed. |
---|
74 | ! |
---|
75 | ! 1960 2016-07-12 16:34:24Z suehring |
---|
76 | ! Scalar surface flux added |
---|
77 | ! |
---|
78 | ! 1849 2016-04-08 11:33:18Z hoffmann |
---|
79 | ! prr moved to arrays_3d |
---|
80 | ! |
---|
81 | ! 1822 2016-04-07 07:49:42Z hoffmann |
---|
82 | ! prr vertical dimensions set to nzb_do to nzt_do. Unused variables deleted. |
---|
83 | ! |
---|
84 | ! 1808 2016-04-05 19:44:00Z raasch |
---|
85 | ! test output removed |
---|
86 | ! |
---|
87 | ! 1783 2016-03-06 18:36:17Z raasch |
---|
88 | ! name change of netcdf routines and module + related changes |
---|
89 | ! |
---|
90 | ! 1745 2016-02-05 13:06:51Z gronemeier |
---|
91 | ! Bugfix: test if time axis limit exceeds moved to point after call of check_open |
---|
92 | ! |
---|
93 | ! 1691 2015-10-26 16:17:44Z maronga |
---|
94 | ! Added output of radiative heating rates for RRTMG |
---|
95 | ! |
---|
96 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
97 | ! Code annotations made doxygen readable |
---|
98 | ! |
---|
99 | ! 1585 2015-04-30 07:05:52Z maronga |
---|
100 | ! Added support for RRTMG |
---|
101 | ! |
---|
102 | ! 1551 2015-03-03 14:18:16Z maronga |
---|
103 | ! Added suppport for land surface model and radiation model output. In the course |
---|
104 | ! of this action, the limits for vertical loops have been changed (from nzb and |
---|
105 | ! nzt+1 to nzb_do and nzt_do, respectively in order to allow soil model output). |
---|
106 | ! Moreover, a new vertical grid zs was introduced. |
---|
107 | ! |
---|
108 | ! 1359 2014-04-11 17:15:14Z hoffmann |
---|
109 | ! New particle structure integrated. |
---|
110 | ! |
---|
111 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
112 | ! REAL constants provided with KIND-attribute |
---|
113 | ! |
---|
114 | ! 1327 2014-03-21 11:00:16Z raasch |
---|
115 | ! parts concerning avs output removed, |
---|
116 | ! -netcdf output queries |
---|
117 | ! |
---|
118 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
119 | ! ONLY-attribute added to USE-statements, |
---|
120 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
121 | ! kinds are defined in new module kinds, |
---|
122 | ! old module precision_kind is removed, |
---|
123 | ! revision history before 2012 removed, |
---|
124 | ! comment fields (!:) to be used for variable explanations added to |
---|
125 | ! all variable declaration statements |
---|
126 | ! |
---|
127 | ! 1318 2014-03-17 13:35:16Z raasch |
---|
128 | ! barrier argument removed from cpu_log, |
---|
129 | ! module interfaces removed |
---|
130 | ! |
---|
131 | ! 1308 2014-03-13 14:58:42Z fricke |
---|
132 | ! Check, if the limit of the time dimension is exceeded for parallel output |
---|
133 | ! To increase the performance for parallel output, the following is done: |
---|
134 | ! - Update of time axis is only done by PE0 |
---|
135 | ! |
---|
136 | ! 1244 2013-10-31 08:16:56Z raasch |
---|
137 | ! Bugfix for index bounds in case of 3d-parallel output |
---|
138 | ! |
---|
139 | ! 1115 2013-03-26 18:16:16Z hoffmann |
---|
140 | ! ql is calculated by calc_liquid_water_content |
---|
141 | ! |
---|
142 | ! 1106 2013-03-04 05:31:38Z raasch |
---|
143 | ! array_kind renamed precision_kind |
---|
144 | ! |
---|
145 | ! 1076 2012-12-05 08:30:18Z hoffmann |
---|
146 | ! Bugfix in output of ql |
---|
147 | ! |
---|
148 | ! 1053 2012-11-13 17:11:03Z hoffmann |
---|
149 | ! +nr, qr, prr, qc and averaged quantities |
---|
150 | ! |
---|
151 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
152 | ! code put under GPL (PALM 3.9) |
---|
153 | ! |
---|
154 | ! 1031 2012-10-19 14:35:30Z raasch |
---|
155 | ! netCDF4 without parallel file support implemented |
---|
156 | ! |
---|
157 | ! 1007 2012-09-19 14:30:36Z franke |
---|
158 | ! Bugfix: missing calculation of ql_vp added |
---|
159 | ! |
---|
160 | ! Revision 1.1 1997/09/03 06:29:36 raasch |
---|
161 | ! Initial revision |
---|
162 | ! |
---|
163 | ! |
---|
164 | ! Description: |
---|
165 | ! ------------ |
---|
166 | !> Output of the 3D-arrays in netCDF and/or AVS format. |
---|
167 | !------------------------------------------------------------------------------! |
---|
168 | SUBROUTINE data_output_3d( av ) |
---|
169 | |
---|
170 | |
---|
171 | USE arrays_3d, & |
---|
172 | ONLY: e, nc, nr, p, pt, prr, q, qc, ql, ql_c, ql_v, qr, rho_ocean, s, & |
---|
173 | sa, tend, u, v, vpt, w |
---|
174 | |
---|
175 | USE averaging |
---|
176 | |
---|
177 | #if defined( __chem ) |
---|
178 | USE chemistry_model_mod, & |
---|
179 | ONLY: chem_data_output_3d |
---|
180 | #endif |
---|
181 | |
---|
182 | USE cloud_parameters, & |
---|
183 | ONLY: l_d_cp, pt_d_t |
---|
184 | |
---|
185 | USE control_parameters, & |
---|
186 | ONLY: air_chemistry, cloud_physics, do3d, do3d_no, do3d_time_count, & |
---|
187 | io_blocks, io_group, land_surface, message_string, & |
---|
188 | ntdim_3d, nz_do3d, & |
---|
189 | psolver, simulated_time, time_since_reference_point, & |
---|
190 | urban_surface, varnamelength |
---|
191 | |
---|
192 | USE cpulog, & |
---|
193 | ONLY: log_point, cpu_log |
---|
194 | |
---|
195 | USE indices, & |
---|
196 | ONLY: nbgp, nx, nxl, nxr, ny, nyn, nys, nzb, wall_flags_0 |
---|
197 | |
---|
198 | USE kinds |
---|
199 | |
---|
200 | USE land_surface_model_mod, & |
---|
201 | ONLY: lsm_data_output_3d, nzb_soil, nzt_soil |
---|
202 | |
---|
203 | #if defined( __netcdf ) |
---|
204 | USE NETCDF |
---|
205 | #endif |
---|
206 | |
---|
207 | USE netcdf_interface, & |
---|
208 | ONLY: fill_value, id_set_3d, id_var_do3d, id_var_time_3d, nc_stat, & |
---|
209 | netcdf_data_format, netcdf_handle_error |
---|
210 | |
---|
211 | USE particle_attributes, & |
---|
212 | ONLY: grid_particles, number_of_particles, particles, & |
---|
213 | particle_advection_start, prt_count |
---|
214 | |
---|
215 | USE pegrid |
---|
216 | |
---|
217 | USE plant_canopy_model_mod, & |
---|
218 | ONLY: pcm_data_output_3d, plant_canopy |
---|
219 | |
---|
220 | USE radiation_model_mod, & |
---|
221 | ONLY: nzub, nzut, radiation, radiation_data_output_3d |
---|
222 | |
---|
223 | USE turbulence_closure_mod, & |
---|
224 | ONLY: tcm_data_output_3d |
---|
225 | |
---|
226 | USE urban_surface_mod, & |
---|
227 | ONLY: usm_data_output_3d |
---|
228 | |
---|
229 | |
---|
230 | IMPLICIT NONE |
---|
231 | |
---|
232 | INTEGER(iwp) :: av !< |
---|
233 | INTEGER(iwp) :: flag_nr !< number of masking flag |
---|
234 | INTEGER(iwp) :: i !< |
---|
235 | INTEGER(iwp) :: if !< |
---|
236 | INTEGER(iwp) :: j !< |
---|
237 | INTEGER(iwp) :: k !< |
---|
238 | INTEGER(iwp) :: n !< |
---|
239 | INTEGER(iwp) :: nzb_do !< vertical lower limit for data output |
---|
240 | INTEGER(iwp) :: nzt_do !< vertical upper limit for data output |
---|
241 | |
---|
242 | LOGICAL :: found !< |
---|
243 | LOGICAL :: resorted !< |
---|
244 | |
---|
245 | REAL(wp) :: mean_r !< |
---|
246 | REAL(wp) :: s_r2 !< |
---|
247 | REAL(wp) :: s_r3 !< |
---|
248 | |
---|
249 | REAL(sp), DIMENSION(:,:,:), ALLOCATABLE :: local_pf !< |
---|
250 | |
---|
251 | REAL(wp), DIMENSION(:,:,:), POINTER :: to_be_resorted !< |
---|
252 | |
---|
253 | CHARACTER (LEN=varnamelength) :: trimvar !< TRIM of output-variable string |
---|
254 | |
---|
255 | ! |
---|
256 | !-- Return, if nothing to output |
---|
257 | IF ( do3d_no(av) == 0 ) RETURN |
---|
258 | |
---|
259 | CALL cpu_log (log_point(14),'data_output_3d','start') |
---|
260 | |
---|
261 | ! |
---|
262 | !-- Open output file. |
---|
263 | !-- For classic or 64bit netCDF output on more than one PE, each PE opens its |
---|
264 | !-- own file and writes the data of its subdomain in binary format. After the |
---|
265 | !-- run, these files are combined to one NetCDF file by combine_plot_fields. |
---|
266 | !-- For netCDF4/HDF5 output, data is written in parallel into one file. |
---|
267 | IF ( netcdf_data_format < 5 ) THEN |
---|
268 | CALL check_open( 30 ) |
---|
269 | IF ( myid == 0 ) CALL check_open( 106+av*10 ) |
---|
270 | ELSE |
---|
271 | CALL check_open( 106+av*10 ) |
---|
272 | ENDIF |
---|
273 | |
---|
274 | ! |
---|
275 | !-- For parallel netcdf output the time axis must be limited. Return, if this |
---|
276 | !-- limit is exceeded. This could be the case, if the simulated time exceeds |
---|
277 | !-- the given end time by the length of the given output interval. |
---|
278 | IF ( netcdf_data_format > 4 ) THEN |
---|
279 | IF ( do3d_time_count(av) + 1 > ntdim_3d(av) ) THEN |
---|
280 | WRITE ( message_string, * ) 'Output of 3d data is not given at t=', & |
---|
281 | simulated_time, '&because the maximum ', & |
---|
282 | 'number of output time levels is ', & |
---|
283 | 'exceeded.' |
---|
284 | CALL message( 'data_output_3d', 'PA0387', 0, 1, 0, 6, 0 ) |
---|
285 | CALL cpu_log( log_point(14), 'data_output_3d', 'stop' ) |
---|
286 | RETURN |
---|
287 | ENDIF |
---|
288 | ENDIF |
---|
289 | |
---|
290 | ! |
---|
291 | !-- Update the netCDF time axis |
---|
292 | !-- In case of parallel output, this is only done by PE0 to increase the |
---|
293 | !-- performance. |
---|
294 | #if defined( __netcdf ) |
---|
295 | do3d_time_count(av) = do3d_time_count(av) + 1 |
---|
296 | IF ( myid == 0 ) THEN |
---|
297 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_time_3d(av), & |
---|
298 | (/ time_since_reference_point /), & |
---|
299 | start = (/ do3d_time_count(av) /), & |
---|
300 | count = (/ 1 /) ) |
---|
301 | CALL netcdf_handle_error( 'data_output_3d', 376 ) |
---|
302 | ENDIF |
---|
303 | #endif |
---|
304 | |
---|
305 | ! |
---|
306 | !-- Loop over all variables to be written. |
---|
307 | if = 1 |
---|
308 | |
---|
309 | DO WHILE ( do3d(av,if)(1:1) /= ' ' ) |
---|
310 | |
---|
311 | ! |
---|
312 | !-- Temporary solution to account for data output within the new urban |
---|
313 | !-- surface model (urban_surface_mod.f90), see also SELECT CASE ( trimvar ). |
---|
314 | !-- Store the array chosen on the temporary array. |
---|
315 | trimvar = TRIM( do3d(av,if) ) |
---|
316 | IF ( urban_surface .AND. trimvar(1:4) == 'usm_' ) THEN |
---|
317 | trimvar = 'usm_output' |
---|
318 | resorted = .TRUE. |
---|
319 | nzb_do = nzub |
---|
320 | nzt_do = nzut |
---|
321 | ELSE |
---|
322 | resorted = .FALSE. |
---|
323 | nzb_do = nzb |
---|
324 | nzt_do = nz_do3d |
---|
325 | ENDIF |
---|
326 | ! |
---|
327 | !-- Set flag to steer output of radiation, land-surface, or user-defined |
---|
328 | !-- quantities |
---|
329 | found = .FALSE. |
---|
330 | ! |
---|
331 | !-- Allocate a temporary array with the desired output dimensions. |
---|
332 | ALLOCATE( local_pf(nxl:nxr,nys:nyn,nzb_do:nzt_do) ) |
---|
333 | ! |
---|
334 | !-- Before each output, set array local_pf to fill value |
---|
335 | local_pf = fill_value |
---|
336 | ! |
---|
337 | !-- Set masking flag for topography for not resorted arrays |
---|
338 | flag_nr = 0 |
---|
339 | |
---|
340 | SELECT CASE ( trimvar ) |
---|
341 | |
---|
342 | CASE ( 'e' ) |
---|
343 | IF ( av == 0 ) THEN |
---|
344 | to_be_resorted => e |
---|
345 | ELSE |
---|
346 | to_be_resorted => e_av |
---|
347 | ENDIF |
---|
348 | |
---|
349 | CASE ( 'lpt' ) |
---|
350 | IF ( av == 0 ) THEN |
---|
351 | to_be_resorted => pt |
---|
352 | ELSE |
---|
353 | to_be_resorted => lpt_av |
---|
354 | ENDIF |
---|
355 | |
---|
356 | CASE ( 'nc' ) |
---|
357 | IF ( av == 0 ) THEN |
---|
358 | to_be_resorted => nc |
---|
359 | ELSE |
---|
360 | to_be_resorted => nc_av |
---|
361 | ENDIF |
---|
362 | |
---|
363 | CASE ( 'nr' ) |
---|
364 | IF ( av == 0 ) THEN |
---|
365 | to_be_resorted => nr |
---|
366 | ELSE |
---|
367 | to_be_resorted => nr_av |
---|
368 | ENDIF |
---|
369 | |
---|
370 | CASE ( 'p' ) |
---|
371 | IF ( av == 0 ) THEN |
---|
372 | IF ( psolver /= 'sor' ) CALL exchange_horiz( p, nbgp ) |
---|
373 | to_be_resorted => p |
---|
374 | ELSE |
---|
375 | IF ( psolver /= 'sor' ) CALL exchange_horiz( p_av, nbgp ) |
---|
376 | to_be_resorted => p_av |
---|
377 | ENDIF |
---|
378 | |
---|
379 | CASE ( 'pc' ) ! particle concentration (requires ghostpoint exchange) |
---|
380 | IF ( av == 0 ) THEN |
---|
381 | IF ( simulated_time >= particle_advection_start ) THEN |
---|
382 | tend = prt_count |
---|
383 | ELSE |
---|
384 | tend = 0.0_wp |
---|
385 | ENDIF |
---|
386 | DO i = nxl, nxr |
---|
387 | DO j = nys, nyn |
---|
388 | DO k = nzb_do, nzt_do |
---|
389 | local_pf(i,j,k) = tend(k,j,i) |
---|
390 | ENDDO |
---|
391 | ENDDO |
---|
392 | ENDDO |
---|
393 | resorted = .TRUE. |
---|
394 | ELSE |
---|
395 | to_be_resorted => pc_av |
---|
396 | ENDIF |
---|
397 | |
---|
398 | CASE ( 'pr' ) ! mean particle radius (effective radius) |
---|
399 | IF ( av == 0 ) THEN |
---|
400 | IF ( simulated_time >= particle_advection_start ) THEN |
---|
401 | DO i = nxl, nxr |
---|
402 | DO j = nys, nyn |
---|
403 | DO k = nzb_do, nzt_do |
---|
404 | number_of_particles = prt_count(k,j,i) |
---|
405 | IF (number_of_particles <= 0) CYCLE |
---|
406 | particles => grid_particles(k,j,i)%particles(1:number_of_particles) |
---|
407 | s_r2 = 0.0_wp |
---|
408 | s_r3 = 0.0_wp |
---|
409 | DO n = 1, number_of_particles |
---|
410 | IF ( particles(n)%particle_mask ) THEN |
---|
411 | s_r2 = s_r2 + particles(n)%radius**2 * & |
---|
412 | particles(n)%weight_factor |
---|
413 | s_r3 = s_r3 + particles(n)%radius**3 * & |
---|
414 | particles(n)%weight_factor |
---|
415 | ENDIF |
---|
416 | ENDDO |
---|
417 | IF ( s_r2 > 0.0_wp ) THEN |
---|
418 | mean_r = s_r3 / s_r2 |
---|
419 | ELSE |
---|
420 | mean_r = 0.0_wp |
---|
421 | ENDIF |
---|
422 | tend(k,j,i) = mean_r |
---|
423 | ENDDO |
---|
424 | ENDDO |
---|
425 | ENDDO |
---|
426 | ELSE |
---|
427 | tend = 0.0_wp |
---|
428 | ENDIF |
---|
429 | DO i = nxl, nxr |
---|
430 | DO j = nys, nyn |
---|
431 | DO k = nzb_do, nzt_do |
---|
432 | local_pf(i,j,k) = tend(k,j,i) |
---|
433 | ENDDO |
---|
434 | ENDDO |
---|
435 | ENDDO |
---|
436 | resorted = .TRUE. |
---|
437 | ELSE |
---|
438 | to_be_resorted => pr_av |
---|
439 | ENDIF |
---|
440 | |
---|
441 | CASE ( 'prr' ) |
---|
442 | IF ( av == 0 ) THEN |
---|
443 | DO i = nxl, nxr |
---|
444 | DO j = nys, nyn |
---|
445 | DO k = nzb_do, nzt_do |
---|
446 | local_pf(i,j,k) = prr(k,j,i) |
---|
447 | ENDDO |
---|
448 | ENDDO |
---|
449 | ENDDO |
---|
450 | ELSE |
---|
451 | DO i = nxl, nxr |
---|
452 | DO j = nys, nyn |
---|
453 | DO k = nzb_do, nzt_do |
---|
454 | local_pf(i,j,k) = prr_av(k,j,i) |
---|
455 | ENDDO |
---|
456 | ENDDO |
---|
457 | ENDDO |
---|
458 | ENDIF |
---|
459 | resorted = .TRUE. |
---|
460 | |
---|
461 | CASE ( 'pt' ) |
---|
462 | IF ( av == 0 ) THEN |
---|
463 | IF ( .NOT. cloud_physics ) THEN |
---|
464 | to_be_resorted => pt |
---|
465 | ELSE |
---|
466 | DO i = nxl, nxr |
---|
467 | DO j = nys, nyn |
---|
468 | DO k = nzb_do, nzt_do |
---|
469 | local_pf(i,j,k) = pt(k,j,i) + l_d_cp * & |
---|
470 | pt_d_t(k) * & |
---|
471 | ql(k,j,i) |
---|
472 | ENDDO |
---|
473 | ENDDO |
---|
474 | ENDDO |
---|
475 | resorted = .TRUE. |
---|
476 | ENDIF |
---|
477 | ELSE |
---|
478 | to_be_resorted => pt_av |
---|
479 | ENDIF |
---|
480 | |
---|
481 | CASE ( 'q' ) |
---|
482 | IF ( av == 0 ) THEN |
---|
483 | to_be_resorted => q |
---|
484 | ELSE |
---|
485 | to_be_resorted => q_av |
---|
486 | ENDIF |
---|
487 | |
---|
488 | CASE ( 'qc' ) |
---|
489 | IF ( av == 0 ) THEN |
---|
490 | to_be_resorted => qc |
---|
491 | ELSE |
---|
492 | to_be_resorted => qc_av |
---|
493 | ENDIF |
---|
494 | |
---|
495 | CASE ( 'ql' ) |
---|
496 | IF ( av == 0 ) THEN |
---|
497 | to_be_resorted => ql |
---|
498 | ELSE |
---|
499 | to_be_resorted => ql_av |
---|
500 | ENDIF |
---|
501 | |
---|
502 | CASE ( 'ql_c' ) |
---|
503 | IF ( av == 0 ) THEN |
---|
504 | to_be_resorted => ql_c |
---|
505 | ELSE |
---|
506 | to_be_resorted => ql_c_av |
---|
507 | ENDIF |
---|
508 | |
---|
509 | CASE ( 'ql_v' ) |
---|
510 | IF ( av == 0 ) THEN |
---|
511 | to_be_resorted => ql_v |
---|
512 | ELSE |
---|
513 | to_be_resorted => ql_v_av |
---|
514 | ENDIF |
---|
515 | |
---|
516 | CASE ( 'ql_vp' ) |
---|
517 | IF ( av == 0 ) THEN |
---|
518 | IF ( simulated_time >= particle_advection_start ) THEN |
---|
519 | DO i = nxl, nxr |
---|
520 | DO j = nys, nyn |
---|
521 | DO k = nzb_do, nzt_do |
---|
522 | number_of_particles = prt_count(k,j,i) |
---|
523 | IF (number_of_particles <= 0) CYCLE |
---|
524 | particles => grid_particles(k,j,i)%particles(1:number_of_particles) |
---|
525 | DO n = 1, number_of_particles |
---|
526 | IF ( particles(n)%particle_mask ) THEN |
---|
527 | tend(k,j,i) = tend(k,j,i) + & |
---|
528 | particles(n)%weight_factor / & |
---|
529 | prt_count(k,j,i) |
---|
530 | ENDIF |
---|
531 | ENDDO |
---|
532 | ENDDO |
---|
533 | ENDDO |
---|
534 | ENDDO |
---|
535 | ELSE |
---|
536 | tend = 0.0_wp |
---|
537 | ENDIF |
---|
538 | DO i = nxl, nxr |
---|
539 | DO j = nys, nyn |
---|
540 | DO k = nzb_do, nzt_do |
---|
541 | local_pf(i,j,k) = tend(k,j,i) |
---|
542 | ENDDO |
---|
543 | ENDDO |
---|
544 | ENDDO |
---|
545 | resorted = .TRUE. |
---|
546 | ELSE |
---|
547 | to_be_resorted => ql_vp_av |
---|
548 | ENDIF |
---|
549 | |
---|
550 | CASE ( 'qr' ) |
---|
551 | IF ( av == 0 ) THEN |
---|
552 | to_be_resorted => qr |
---|
553 | ELSE |
---|
554 | to_be_resorted => qr_av |
---|
555 | ENDIF |
---|
556 | |
---|
557 | CASE ( 'qv' ) |
---|
558 | IF ( av == 0 ) THEN |
---|
559 | DO i = nxl, nxr |
---|
560 | DO j = nys, nyn |
---|
561 | DO k = nzb_do, nzt_do |
---|
562 | local_pf(i,j,k) = q(k,j,i) - ql(k,j,i) |
---|
563 | ENDDO |
---|
564 | ENDDO |
---|
565 | ENDDO |
---|
566 | resorted = .TRUE. |
---|
567 | ELSE |
---|
568 | to_be_resorted => qv_av |
---|
569 | ENDIF |
---|
570 | |
---|
571 | CASE ( 'rho_ocean' ) |
---|
572 | IF ( av == 0 ) THEN |
---|
573 | to_be_resorted => rho_ocean |
---|
574 | ELSE |
---|
575 | to_be_resorted => rho_ocean_av |
---|
576 | ENDIF |
---|
577 | |
---|
578 | CASE ( 's' ) |
---|
579 | IF ( av == 0 ) THEN |
---|
580 | to_be_resorted => s |
---|
581 | ELSE |
---|
582 | to_be_resorted => s_av |
---|
583 | ENDIF |
---|
584 | |
---|
585 | CASE ( 'sa' ) |
---|
586 | IF ( av == 0 ) THEN |
---|
587 | to_be_resorted => sa |
---|
588 | ELSE |
---|
589 | to_be_resorted => sa_av |
---|
590 | ENDIF |
---|
591 | |
---|
592 | CASE ( 'u' ) |
---|
593 | flag_nr = 1 |
---|
594 | IF ( av == 0 ) THEN |
---|
595 | to_be_resorted => u |
---|
596 | ELSE |
---|
597 | to_be_resorted => u_av |
---|
598 | ENDIF |
---|
599 | |
---|
600 | CASE ( 'v' ) |
---|
601 | flag_nr = 2 |
---|
602 | IF ( av == 0 ) THEN |
---|
603 | to_be_resorted => v |
---|
604 | ELSE |
---|
605 | to_be_resorted => v_av |
---|
606 | ENDIF |
---|
607 | |
---|
608 | CASE ( 'vpt' ) |
---|
609 | IF ( av == 0 ) THEN |
---|
610 | to_be_resorted => vpt |
---|
611 | ELSE |
---|
612 | to_be_resorted => vpt_av |
---|
613 | ENDIF |
---|
614 | |
---|
615 | CASE ( 'w' ) |
---|
616 | flag_nr = 3 |
---|
617 | IF ( av == 0 ) THEN |
---|
618 | to_be_resorted => w |
---|
619 | ELSE |
---|
620 | to_be_resorted => w_av |
---|
621 | ENDIF |
---|
622 | ! |
---|
623 | !-- Block of urban surface model outputs |
---|
624 | CASE ( 'usm_output' ) |
---|
625 | CALL usm_data_output_3d( av, do3d(av,if), found, local_pf, & |
---|
626 | nzb_do, nzt_do ) |
---|
627 | |
---|
628 | CASE DEFAULT |
---|
629 | |
---|
630 | ! |
---|
631 | !-- Land surface quantity |
---|
632 | IF ( land_surface ) THEN |
---|
633 | ! |
---|
634 | !-- For soil model quantities, it is required to re-allocate local_pf |
---|
635 | nzb_do = nzb_soil |
---|
636 | nzt_do = nzt_soil |
---|
637 | |
---|
638 | DEALLOCATE ( local_pf ) |
---|
639 | ALLOCATE( local_pf(nxl:nxr,nys:nyn,nzb_do:nzt_do) ) |
---|
640 | local_pf = fill_value |
---|
641 | |
---|
642 | CALL lsm_data_output_3d( av, do3d(av,if), found, local_pf ) |
---|
643 | resorted = .TRUE. |
---|
644 | |
---|
645 | ! |
---|
646 | !-- If no soil model variable was found, re-allocate local_pf |
---|
647 | IF ( .NOT. found ) THEN |
---|
648 | nzb_do = nzb |
---|
649 | nzt_do = nz_do3d |
---|
650 | |
---|
651 | DEALLOCATE ( local_pf ) |
---|
652 | ALLOCATE( local_pf(nxl:nxr,nys:nyn,nzb_do:nzt_do) ) |
---|
653 | ENDIF |
---|
654 | |
---|
655 | ENDIF |
---|
656 | |
---|
657 | IF ( .NOT. found ) THEN |
---|
658 | CALL tcm_data_output_3d( av, do3d(av,if), found, local_pf ) |
---|
659 | resorted = .TRUE. |
---|
660 | ENDIF |
---|
661 | |
---|
662 | ! |
---|
663 | !-- Radiation quantity |
---|
664 | IF ( .NOT. found .AND. radiation ) THEN |
---|
665 | CALL radiation_data_output_3d( av, do3d(av,if), found, & |
---|
666 | local_pf ) |
---|
667 | resorted = .TRUE. |
---|
668 | ENDIF |
---|
669 | |
---|
670 | ! |
---|
671 | !-- Chemistry model output |
---|
672 | #if defined( __chem ) |
---|
673 | IF ( .NOT. found .AND. air_chemistry ) THEN |
---|
674 | CALL chem_data_output_3d( av, do3d(av,if), found, & |
---|
675 | local_pf ) |
---|
676 | resorted = .TRUE. |
---|
677 | ENDIF |
---|
678 | #endif |
---|
679 | |
---|
680 | ! |
---|
681 | !-- Plant canopy model output |
---|
682 | IF ( .NOT. found .AND. plant_canopy ) THEN |
---|
683 | CALL pcm_data_output_3d( av, do3d(av,if), found, local_pf, & |
---|
684 | fill_value ) |
---|
685 | resorted = .TRUE. |
---|
686 | ENDIF |
---|
687 | |
---|
688 | ! |
---|
689 | !-- User defined quantity |
---|
690 | IF ( .NOT. found ) THEN |
---|
691 | CALL user_data_output_3d( av, do3d(av,if), found, local_pf, & |
---|
692 | nzb_do, nzt_do ) |
---|
693 | resorted = .TRUE. |
---|
694 | ENDIF |
---|
695 | |
---|
696 | IF ( .NOT. found ) THEN |
---|
697 | message_string = 'no output available for: ' // & |
---|
698 | TRIM( do3d(av,if) ) |
---|
699 | CALL message( 'data_output_3d', 'PA0182', 0, 0, 0, 6, 0 ) |
---|
700 | ENDIF |
---|
701 | |
---|
702 | END SELECT |
---|
703 | |
---|
704 | ! |
---|
705 | !-- Resort the array to be output, if not done above |
---|
706 | IF ( .NOT. resorted ) THEN |
---|
707 | DO i = nxl, nxr |
---|
708 | DO j = nys, nyn |
---|
709 | DO k = nzb_do, nzt_do |
---|
710 | local_pf(i,j,k) = MERGE( & |
---|
711 | to_be_resorted(k,j,i), & |
---|
712 | REAL( fill_value, KIND = wp ), & |
---|
713 | BTEST( wall_flags_0(k,j,i), flag_nr ) ) |
---|
714 | ENDDO |
---|
715 | ENDDO |
---|
716 | ENDDO |
---|
717 | ENDIF |
---|
718 | |
---|
719 | ! |
---|
720 | !-- Output of the 3D-array |
---|
721 | #if defined( __parallel ) |
---|
722 | IF ( netcdf_data_format < 5 ) THEN |
---|
723 | ! |
---|
724 | !-- Non-parallel netCDF output. Data is output in parallel in |
---|
725 | !-- FORTRAN binary format here, and later collected into one file by |
---|
726 | !-- combine_plot_fields |
---|
727 | IF ( myid == 0 ) THEN |
---|
728 | WRITE ( 30 ) time_since_reference_point, & |
---|
729 | do3d_time_count(av), av |
---|
730 | ENDIF |
---|
731 | DO i = 0, io_blocks-1 |
---|
732 | IF ( i == io_group ) THEN |
---|
733 | WRITE ( 30 ) nxl, nxr, nys, nyn, nzb_do, nzt_do |
---|
734 | WRITE ( 30 ) local_pf(:,:,nzb_do:nzt_do) |
---|
735 | ENDIF |
---|
736 | |
---|
737 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
738 | |
---|
739 | ENDDO |
---|
740 | |
---|
741 | ELSE |
---|
742 | #if defined( __netcdf ) |
---|
743 | ! |
---|
744 | !-- Parallel output in netCDF4/HDF5 format. |
---|
745 | ! IF ( nxr == nx .AND. nyn /= ny ) THEN |
---|
746 | ! nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
747 | ! local_pf(nxl:nxr+1,nys:nyn,nzb_do:nzt_do), & |
---|
748 | ! start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
749 | ! count = (/ nxr-nxl+2, nyn-nys+1, nzt_do-nzb_do+1, 1 /) ) |
---|
750 | ! ELSEIF ( nxr /= nx .AND. nyn == ny ) THEN |
---|
751 | ! nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
752 | ! local_pf(nxl:nxr,nys:nyn+1,nzb_do:nzt_do), & |
---|
753 | ! start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
754 | ! count = (/ nxr-nxl+1, nyn-nys+2, nzt_do-nzb_do+1, 1 /) ) |
---|
755 | ! ELSEIF ( nxr == nx .AND. nyn == ny ) THEN |
---|
756 | ! nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
757 | ! local_pf(nxl:nxr+1,nys:nyn+1,nzb_do:nzt_do ), & |
---|
758 | ! start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
759 | ! count = (/ nxr-nxl+2, nyn-nys+2, nzt_do-nzb_do+1, 1 /) ) |
---|
760 | ! ELSE |
---|
761 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
762 | local_pf(nxl:nxr,nys:nyn,nzb_do:nzt_do), & |
---|
763 | start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
764 | count = (/ nxr-nxl+1, nyn-nys+1, nzt_do-nzb_do+1, 1 /) ) |
---|
765 | ! ENDIF |
---|
766 | CALL netcdf_handle_error( 'data_output_3d', 386 ) |
---|
767 | #endif |
---|
768 | ENDIF |
---|
769 | #else |
---|
770 | #if defined( __netcdf ) |
---|
771 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
772 | local_pf(nxl:nxr+1,nys:nyn+1,nzb_do:nzt_do), & |
---|
773 | start = (/ 1, 1, 1, do3d_time_count(av) /), & |
---|
774 | count = (/ nx+1, ny+1, nzt_do-nzb_do+1, 1 /) ) |
---|
775 | CALL netcdf_handle_error( 'data_output_3d', 446 ) |
---|
776 | #endif |
---|
777 | #endif |
---|
778 | |
---|
779 | if = if + 1 |
---|
780 | |
---|
781 | ! |
---|
782 | !-- Deallocate temporary array |
---|
783 | DEALLOCATE ( local_pf ) |
---|
784 | |
---|
785 | ENDDO |
---|
786 | |
---|
787 | CALL cpu_log( log_point(14), 'data_output_3d', 'stop' ) |
---|
788 | |
---|
789 | ! |
---|
790 | !-- Formats. |
---|
791 | 3300 FORMAT ('variable ',I4,' file=',A,' filetype=unformatted skip=',I12/ & |
---|
792 | 'label = ',A,A) |
---|
793 | |
---|
794 | END SUBROUTINE data_output_3d |
---|