1 | SUBROUTINE header |
---|
2 | |
---|
3 | !--------------------------------------------------------------------------------! |
---|
4 | ! This file is part of PALM. |
---|
5 | ! |
---|
6 | ! PALM is free software: you can redistribute it and/or modify it under the terms |
---|
7 | ! of the GNU General Public License as published by the Free Software Foundation, |
---|
8 | ! either version 3 of the License, or (at your option) any later 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-2014 Leibniz Universitaet Hannover |
---|
18 | !--------------------------------------------------------------------------------! |
---|
19 | ! |
---|
20 | ! Current revisions: |
---|
21 | ! ----------------- |
---|
22 | ! |
---|
23 | ! |
---|
24 | ! Former revisions: |
---|
25 | ! ----------------- |
---|
26 | ! $Id: header.f90 1576 2015-03-27 10:23:30Z suehring $ |
---|
27 | ! |
---|
28 | ! 1575 2015-03-27 09:56:27Z raasch |
---|
29 | ! adjustments for psolver-queries, output of seed_follows_topography |
---|
30 | ! |
---|
31 | ! 1560 2015-03-06 10:48:54Z keck |
---|
32 | ! output for recycling y shift |
---|
33 | ! |
---|
34 | ! 1557 2015-03-05 16:43:04Z suehring |
---|
35 | ! output for monotonic limiter |
---|
36 | ! |
---|
37 | ! 1551 2015-03-03 14:18:16Z maronga |
---|
38 | ! Added informal output for land surface model and radiation model. Removed typo. |
---|
39 | ! |
---|
40 | ! 1496 2014-12-02 17:25:50Z maronga |
---|
41 | ! Renamed: "radiation -> "cloud_top_radiation" |
---|
42 | ! |
---|
43 | ! 1484 2014-10-21 10:53:05Z kanani |
---|
44 | ! Changes due to new module structure of the plant canopy model: |
---|
45 | ! module plant_canopy_model_mod and output for new canopy model parameters |
---|
46 | ! (alpha_lad, beta_lad, lai_beta,...) added, |
---|
47 | ! drag_coefficient, leaf_surface_concentration and scalar_exchange_coefficient |
---|
48 | ! renamed to canopy_drag_coeff, leaf_surface_conc and leaf_scalar_exch_coeff, |
---|
49 | ! learde renamed leaf_area_density. |
---|
50 | ! Bugfix: DO-WHILE-loop for lad header information additionally restricted |
---|
51 | ! by maximum number of gradient levels (currently 10) |
---|
52 | ! |
---|
53 | ! 1482 2014-10-18 12:34:45Z raasch |
---|
54 | ! information about calculated or predefined virtual processor topology adjusted |
---|
55 | ! |
---|
56 | ! 1468 2014-09-24 14:06:57Z maronga |
---|
57 | ! Adapted for use on up to 6-digit processor cores |
---|
58 | ! |
---|
59 | ! 1429 2014-07-15 12:53:45Z knoop |
---|
60 | ! header exended to provide ensemble_member_nr if specified |
---|
61 | ! |
---|
62 | ! 1376 2014-04-26 11:21:22Z boeske |
---|
63 | ! Correction of typos |
---|
64 | ! |
---|
65 | ! 1365 2014-04-22 15:03:56Z boeske |
---|
66 | ! New section 'Large scale forcing and nudging': |
---|
67 | ! output of large scale forcing and nudging information, |
---|
68 | ! new section for initial profiles created |
---|
69 | ! |
---|
70 | ! 1359 2014-04-11 17:15:14Z hoffmann |
---|
71 | ! dt_sort_particles removed |
---|
72 | ! |
---|
73 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
74 | ! REAL constants provided with KIND-attribute |
---|
75 | ! |
---|
76 | ! 1327 2014-03-21 11:00:16Z raasch |
---|
77 | ! parts concerning iso2d and avs output removed, |
---|
78 | ! -netcdf output queries |
---|
79 | ! |
---|
80 | ! 1324 2014-03-21 09:13:16Z suehring |
---|
81 | ! Bugfix: module spectrum added |
---|
82 | ! |
---|
83 | ! 1322 2014-03-20 16:38:49Z raasch |
---|
84 | ! REAL functions provided with KIND-attribute, |
---|
85 | ! some REAL constants defined as wp-kind |
---|
86 | ! |
---|
87 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
88 | ! ONLY-attribute added to USE-statements, |
---|
89 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
90 | ! kinds are defined in new module kinds, |
---|
91 | ! revision history before 2012 removed, |
---|
92 | ! comment fields (!:) to be used for variable explanations added to |
---|
93 | ! all variable declaration statements |
---|
94 | ! |
---|
95 | ! 1308 2014-03-13 14:58:42Z fricke |
---|
96 | ! output of the fixed number of output time levels |
---|
97 | ! output_format adjusted for masked data if netcdf_data_format > 5 |
---|
98 | ! |
---|
99 | ! 1299 2014-03-06 13:15:21Z heinze |
---|
100 | ! output for using large_scale subsidence in combination |
---|
101 | ! with large_scale_forcing |
---|
102 | ! reformatting, more detailed explanations |
---|
103 | ! |
---|
104 | ! 1241 2013-10-30 11:36:58Z heinze |
---|
105 | ! output for nudging + large scale forcing from external file |
---|
106 | ! |
---|
107 | ! 1216 2013-08-26 09:31:42Z raasch |
---|
108 | ! output for transpose_compute_overlap |
---|
109 | ! |
---|
110 | ! 1212 2013-08-15 08:46:27Z raasch |
---|
111 | ! output for poisfft_hybrid removed |
---|
112 | ! |
---|
113 | ! 1179 2013-06-14 05:57:58Z raasch |
---|
114 | ! output of reference_state, use_reference renamed use_single_reference_value |
---|
115 | ! |
---|
116 | ! 1159 2013-05-21 11:58:22Z fricke |
---|
117 | ! +use_cmax |
---|
118 | ! |
---|
119 | ! 1115 2013-03-26 18:16:16Z hoffmann |
---|
120 | ! descriptions for Seifert-Beheng-cloud-physics-scheme added |
---|
121 | ! |
---|
122 | ! 1111 2013-03-08 23:54:10Z raasch |
---|
123 | ! output of accelerator board information |
---|
124 | ! ibc_p_b = 2 removed |
---|
125 | ! |
---|
126 | ! 1108 2013-03-05 07:03:32Z raasch |
---|
127 | ! bugfix for r1106 |
---|
128 | ! |
---|
129 | ! 1106 2013-03-04 05:31:38Z raasch |
---|
130 | ! some format changes for coupled runs |
---|
131 | ! |
---|
132 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
133 | ! unused variables removed |
---|
134 | ! |
---|
135 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
136 | ! code put under GPL (PALM 3.9) |
---|
137 | ! |
---|
138 | ! 1031 2012-10-19 14:35:30Z raasch |
---|
139 | ! output of netCDF data format modified |
---|
140 | ! |
---|
141 | ! 1015 2012-09-27 09:23:24Z raasch |
---|
142 | ! output of Adjustment of mixing length to the Prandtl mixing length at first |
---|
143 | ! grid point above ground removed |
---|
144 | ! |
---|
145 | ! 1003 2012-09-14 14:35:53Z raasch |
---|
146 | ! output of information about equal/unequal subdomain size removed |
---|
147 | ! |
---|
148 | ! 1001 2012-09-13 14:08:46Z raasch |
---|
149 | ! all actions concerning leapfrog- and upstream-spline-scheme removed |
---|
150 | ! |
---|
151 | ! 978 2012-08-09 08:28:32Z fricke |
---|
152 | ! -km_damp_max, outflow_damping_width |
---|
153 | ! +pt_damping_factor, pt_damping_width |
---|
154 | ! +z0h |
---|
155 | ! |
---|
156 | ! 964 2012-07-26 09:14:24Z raasch |
---|
157 | ! output of profil-related quantities removed |
---|
158 | ! |
---|
159 | ! 940 2012-07-09 14:31:00Z raasch |
---|
160 | ! Output in case of simulations for pure neutral stratification (no pt-equation |
---|
161 | ! solved) |
---|
162 | ! |
---|
163 | ! 927 2012-06-06 19:15:04Z raasch |
---|
164 | ! output of masking_method for mg-solver |
---|
165 | ! |
---|
166 | ! 868 2012-03-28 12:21:07Z raasch |
---|
167 | ! translation velocity in Galilean transformation changed to 0.6 * ug |
---|
168 | ! |
---|
169 | ! 833 2012-02-22 08:55:55Z maronga |
---|
170 | ! Adjusted format for leaf area density |
---|
171 | ! |
---|
172 | ! 828 2012-02-21 12:00:36Z raasch |
---|
173 | ! output of dissipation_classes + radius_classes |
---|
174 | ! |
---|
175 | ! 825 2012-02-19 03:03:44Z raasch |
---|
176 | ! Output of cloud physics parameters/quantities complemented and restructured |
---|
177 | ! |
---|
178 | ! Revision 1.1 1997/08/11 06:17:20 raasch |
---|
179 | ! Initial revision |
---|
180 | ! |
---|
181 | ! |
---|
182 | ! Description: |
---|
183 | ! ------------ |
---|
184 | ! Writing a header with all important information about the actual run. |
---|
185 | ! This subroutine is called three times, two times at the beginning |
---|
186 | ! (writing information on files RUN_CONTROL and HEADER) and one time at the |
---|
187 | ! end of the run, then writing additional information about CPU-usage on file |
---|
188 | ! header. |
---|
189 | !-----------------------------------------------------------------------------! |
---|
190 | |
---|
191 | USE arrays_3d, & |
---|
192 | ONLY: pt_init, qsws, q_init, sa_init, shf, ug, vg, w_subs, zu |
---|
193 | |
---|
194 | USE control_parameters |
---|
195 | |
---|
196 | USE cloud_parameters, & |
---|
197 | ONLY: cp, curvature_solution_effects, c_sedimentation, & |
---|
198 | limiter_sedimentation, l_v, nc_const, r_d, ventilation_effect |
---|
199 | |
---|
200 | USE cpulog, & |
---|
201 | ONLY: log_point_s |
---|
202 | |
---|
203 | USE dvrp_variables, & |
---|
204 | ONLY: use_seperate_pe_for_dvrp_output |
---|
205 | |
---|
206 | USE grid_variables, & |
---|
207 | ONLY: dx, dy |
---|
208 | |
---|
209 | USE indices, & |
---|
210 | ONLY: mg_loc_ind, nnx, nny, nnz, nx, ny, nxl_mg, nxr_mg, nyn_mg, & |
---|
211 | nys_mg, nzt, nzt_mg |
---|
212 | |
---|
213 | USE kinds |
---|
214 | |
---|
215 | USE land_surface_model_mod, & |
---|
216 | ONLY: conserve_water_content, dewfall, land_surface, nzb_soil, & |
---|
217 | nzt_soil, root_fraction, soil_moisture, soil_temperature, & |
---|
218 | soil_type, soil_type_name, veg_type, veg_type_name, zs |
---|
219 | |
---|
220 | USE model_1d, & |
---|
221 | ONLY: damp_level_ind_1d, dt_pr_1d, dt_run_control_1d, end_time_1d |
---|
222 | |
---|
223 | USE particle_attributes, & |
---|
224 | ONLY: bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, collision_kernel, & |
---|
225 | density_ratio, dissipation_classes, dt_min_part, dt_prel, & |
---|
226 | dt_write_particle_data, end_time_prel, & |
---|
227 | maximum_number_of_tailpoints, maximum_tailpoint_age, & |
---|
228 | minimum_tailpoint_distance, number_of_particle_groups, & |
---|
229 | particle_advection, particle_advection_start, & |
---|
230 | particles_per_point, pdx, pdy, pdz, psb, psl, psn, psr, pss, & |
---|
231 | pst, radius, radius_classes, random_start_position, & |
---|
232 | seed_follows_topography, & |
---|
233 | total_number_of_particles, use_particle_tails, & |
---|
234 | use_sgs_for_particles, total_number_of_tails, & |
---|
235 | vertical_particle_advection, write_particle_statistics |
---|
236 | |
---|
237 | USE pegrid |
---|
238 | |
---|
239 | USE plant_canopy_model_mod, & |
---|
240 | ONLY: alpha_lad, beta_lad, calc_beta_lad_profile, canopy_drag_coeff, & |
---|
241 | canopy_mode, cthf, lad, lad_surface, lad_vertical_gradient, & |
---|
242 | lad_vertical_gradient_level, lad_vertical_gradient_level_ind, & |
---|
243 | lai_beta, leaf_scalar_exch_coeff, leaf_surface_conc, pch_index, & |
---|
244 | plant_canopy |
---|
245 | |
---|
246 | USE radiation_model_mod, & |
---|
247 | ONLY: albedo, day_init, dt_radiation, lambda, net_radiation, & |
---|
248 | radiation, radiation_scheme, time_utc_init |
---|
249 | |
---|
250 | USE spectrum, & |
---|
251 | ONLY: comp_spectra_level, data_output_sp, plot_spectra_level, & |
---|
252 | spectra_direction |
---|
253 | |
---|
254 | IMPLICIT NONE |
---|
255 | |
---|
256 | CHARACTER (LEN=1) :: prec !: |
---|
257 | |
---|
258 | CHARACTER (LEN=2) :: do2d_mode !: |
---|
259 | |
---|
260 | CHARACTER (LEN=5) :: section_chr !: |
---|
261 | |
---|
262 | CHARACTER (LEN=10) :: coor_chr !: |
---|
263 | CHARACTER (LEN=10) :: host_chr !: |
---|
264 | |
---|
265 | CHARACTER (LEN=16) :: begin_chr !: |
---|
266 | |
---|
267 | CHARACTER (LEN=26) :: ver_rev !: |
---|
268 | |
---|
269 | CHARACTER (LEN=40) :: output_format !: |
---|
270 | |
---|
271 | CHARACTER (LEN=70) :: char1 !: |
---|
272 | CHARACTER (LEN=70) :: char2 !: |
---|
273 | CHARACTER (LEN=70) :: dopr_chr !: |
---|
274 | CHARACTER (LEN=70) :: do2d_xy !: |
---|
275 | CHARACTER (LEN=70) :: do2d_xz !: |
---|
276 | CHARACTER (LEN=70) :: do2d_yz !: |
---|
277 | CHARACTER (LEN=70) :: do3d_chr !: |
---|
278 | CHARACTER (LEN=70) :: domask_chr !: |
---|
279 | CHARACTER (LEN=70) :: run_classification !: |
---|
280 | |
---|
281 | CHARACTER (LEN=85) :: roben !: |
---|
282 | CHARACTER (LEN=85) :: runten !: |
---|
283 | |
---|
284 | CHARACTER (LEN=86) :: coordinates !: |
---|
285 | CHARACTER (LEN=86) :: gradients !: |
---|
286 | CHARACTER (LEN=86) :: leaf_area_density !: |
---|
287 | CHARACTER (LEN=86) :: roots !: |
---|
288 | CHARACTER (LEN=86) :: slices !: |
---|
289 | CHARACTER (LEN=86) :: temperatures !: |
---|
290 | CHARACTER (LEN=86) :: ugcomponent !: |
---|
291 | CHARACTER (LEN=86) :: vgcomponent !: |
---|
292 | |
---|
293 | CHARACTER (LEN=1), DIMENSION(1:3) :: dir = (/ 'x', 'y', 'z' /) !: |
---|
294 | |
---|
295 | INTEGER(iwp) :: av !: |
---|
296 | INTEGER(iwp) :: bh !: |
---|
297 | INTEGER(iwp) :: blx !: |
---|
298 | INTEGER(iwp) :: bly !: |
---|
299 | INTEGER(iwp) :: bxl !: |
---|
300 | INTEGER(iwp) :: bxr !: |
---|
301 | INTEGER(iwp) :: byn !: |
---|
302 | INTEGER(iwp) :: bys !: |
---|
303 | INTEGER(iwp) :: ch !: |
---|
304 | INTEGER(iwp) :: count !: |
---|
305 | INTEGER(iwp) :: cwx !: |
---|
306 | INTEGER(iwp) :: cwy !: |
---|
307 | INTEGER(iwp) :: cxl !: |
---|
308 | INTEGER(iwp) :: cxr !: |
---|
309 | INTEGER(iwp) :: cyn !: |
---|
310 | INTEGER(iwp) :: cys !: |
---|
311 | INTEGER(iwp) :: dim !: |
---|
312 | INTEGER(iwp) :: i !: |
---|
313 | INTEGER(iwp) :: io !: |
---|
314 | INTEGER(iwp) :: j !: |
---|
315 | INTEGER(iwp) :: k !: |
---|
316 | INTEGER(iwp) :: l !: |
---|
317 | INTEGER(iwp) :: ll !: |
---|
318 | INTEGER(iwp) :: mpi_type !: |
---|
319 | |
---|
320 | REAL(wp) :: canopy_height !: canopy height (in m) |
---|
321 | REAL(wp) :: cpuseconds_per_simulated_second !: |
---|
322 | |
---|
323 | ! |
---|
324 | !-- Open the output file. At the end of the simulation, output is directed |
---|
325 | !-- to unit 19. |
---|
326 | IF ( ( runnr == 0 .OR. force_print_header ) .AND. & |
---|
327 | .NOT. simulated_time_at_begin /= simulated_time ) THEN |
---|
328 | io = 15 ! header output on file RUN_CONTROL |
---|
329 | ELSE |
---|
330 | io = 19 ! header output on file HEADER |
---|
331 | ENDIF |
---|
332 | CALL check_open( io ) |
---|
333 | |
---|
334 | ! |
---|
335 | !-- At the end of the run, output file (HEADER) will be rewritten with |
---|
336 | !-- new information |
---|
337 | IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 ) |
---|
338 | |
---|
339 | ! |
---|
340 | !-- Determine kind of model run |
---|
341 | IF ( TRIM( initializing_actions ) == 'read_restart_data' ) THEN |
---|
342 | run_classification = '3D - restart run' |
---|
343 | ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' ) THEN |
---|
344 | run_classification = '3D - run with cyclic fill of 3D - prerun data' |
---|
345 | ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 ) THEN |
---|
346 | run_classification = '3D - run without 1D - prerun' |
---|
347 | ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 ) THEN |
---|
348 | run_classification = '3D - run with 1D - prerun' |
---|
349 | ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 ) THEN |
---|
350 | run_classification = '3D - run initialized by user' |
---|
351 | ELSE |
---|
352 | message_string = ' unknown action(s): ' // TRIM( initializing_actions ) |
---|
353 | CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 ) |
---|
354 | ENDIF |
---|
355 | IF ( ocean ) THEN |
---|
356 | run_classification = 'ocean - ' // run_classification |
---|
357 | ELSE |
---|
358 | run_classification = 'atmosphere - ' // run_classification |
---|
359 | ENDIF |
---|
360 | |
---|
361 | ! |
---|
362 | !-- Run-identification, date, time, host |
---|
363 | host_chr = host(1:10) |
---|
364 | ver_rev = TRIM( version ) // ' ' // TRIM( revision ) |
---|
365 | WRITE ( io, 100 ) ver_rev, TRIM( run_classification ) |
---|
366 | IF ( TRIM( coupling_mode ) /= 'uncoupled' ) THEN |
---|
367 | #if defined( __mpi2 ) |
---|
368 | mpi_type = 2 |
---|
369 | #else |
---|
370 | mpi_type = 1 |
---|
371 | #endif |
---|
372 | WRITE ( io, 101 ) mpi_type, coupling_mode |
---|
373 | ENDIF |
---|
374 | #if defined( __parallel ) |
---|
375 | IF ( coupling_start_time /= 0.0_wp ) THEN |
---|
376 | IF ( coupling_start_time > simulated_time_at_begin ) THEN |
---|
377 | WRITE ( io, 109 ) |
---|
378 | ELSE |
---|
379 | WRITE ( io, 114 ) |
---|
380 | ENDIF |
---|
381 | ENDIF |
---|
382 | #endif |
---|
383 | IF ( ensemble_member_nr /= 0 ) THEN |
---|
384 | WRITE ( io, 512 ) run_date, run_identifier, run_time, runnr, & |
---|
385 | ADJUSTR( host_chr ), ensemble_member_nr |
---|
386 | ELSE |
---|
387 | WRITE ( io, 102 ) run_date, run_identifier, run_time, runnr, & |
---|
388 | ADJUSTR( host_chr ) |
---|
389 | ENDIF |
---|
390 | #if defined( __parallel ) |
---|
391 | IF ( npex == -1 .AND. npey == -1 ) THEN |
---|
392 | char1 = 'calculated' |
---|
393 | ELSE |
---|
394 | char1 = 'predefined' |
---|
395 | ENDIF |
---|
396 | IF ( threads_per_task == 1 ) THEN |
---|
397 | WRITE ( io, 103 ) numprocs, pdims(1), pdims(2), TRIM( char1 ) |
---|
398 | ELSE |
---|
399 | WRITE ( io, 104 ) numprocs*threads_per_task, numprocs, & |
---|
400 | threads_per_task, pdims(1), pdims(2), TRIM( char1 ) |
---|
401 | ENDIF |
---|
402 | IF ( num_acc_per_node /= 0 ) WRITE ( io, 117 ) num_acc_per_node |
---|
403 | IF ( ( host(1:3) == 'ibm' .OR. host(1:3) == 'nec' .OR. & |
---|
404 | host(1:2) == 'lc' .OR. host(1:3) == 'dec' ) .AND. & |
---|
405 | npex == -1 .AND. pdims(2) == 1 ) & |
---|
406 | THEN |
---|
407 | WRITE ( io, 106 ) |
---|
408 | ELSEIF ( pdims(2) == 1 ) THEN |
---|
409 | WRITE ( io, 107 ) 'x' |
---|
410 | ELSEIF ( pdims(1) == 1 ) THEN |
---|
411 | WRITE ( io, 107 ) 'y' |
---|
412 | ENDIF |
---|
413 | IF ( use_seperate_pe_for_dvrp_output ) WRITE ( io, 105 ) |
---|
414 | IF ( numprocs /= maximum_parallel_io_streams ) THEN |
---|
415 | WRITE ( io, 108 ) maximum_parallel_io_streams |
---|
416 | ENDIF |
---|
417 | #else |
---|
418 | IF ( num_acc_per_node /= 0 ) WRITE ( io, 120 ) num_acc_per_node |
---|
419 | #endif |
---|
420 | WRITE ( io, 99 ) |
---|
421 | |
---|
422 | ! |
---|
423 | !-- Numerical schemes |
---|
424 | WRITE ( io, 110 ) |
---|
425 | IF ( psolver(1:7) == 'poisfft' ) THEN |
---|
426 | WRITE ( io, 111 ) TRIM( fft_method ) |
---|
427 | IF ( transpose_compute_overlap ) WRITE( io, 115 ) |
---|
428 | ELSEIF ( psolver == 'sor' ) THEN |
---|
429 | WRITE ( io, 112 ) nsor_ini, nsor, omega_sor |
---|
430 | ELSEIF ( psolver(1:9) == 'multigrid' ) THEN |
---|
431 | WRITE ( io, 135 ) TRIM(psolver), cycle_mg, maximum_grid_level, ngsrb |
---|
432 | IF ( mg_cycles == -1 ) THEN |
---|
433 | WRITE ( io, 140 ) residual_limit |
---|
434 | ELSE |
---|
435 | WRITE ( io, 141 ) mg_cycles |
---|
436 | ENDIF |
---|
437 | IF ( mg_switch_to_pe0_level == 0 ) THEN |
---|
438 | WRITE ( io, 136 ) nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, & |
---|
439 | nzt_mg(1) |
---|
440 | ELSEIF ( mg_switch_to_pe0_level /= -1 ) THEN |
---|
441 | WRITE ( io, 137 ) mg_switch_to_pe0_level, & |
---|
442 | mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, & |
---|
443 | mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, & |
---|
444 | nzt_mg(mg_switch_to_pe0_level), & |
---|
445 | nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, & |
---|
446 | nzt_mg(1) |
---|
447 | ENDIF |
---|
448 | IF ( masking_method ) WRITE ( io, 144 ) |
---|
449 | ENDIF |
---|
450 | IF ( call_psolver_at_all_substeps .AND. timestep_scheme(1:5) == 'runge' ) & |
---|
451 | THEN |
---|
452 | WRITE ( io, 142 ) |
---|
453 | ENDIF |
---|
454 | |
---|
455 | IF ( momentum_advec == 'pw-scheme' ) THEN |
---|
456 | WRITE ( io, 113 ) |
---|
457 | ELSEIF (momentum_advec == 'ws-scheme' ) THEN |
---|
458 | WRITE ( io, 503 ) |
---|
459 | ENDIF |
---|
460 | IF ( scalar_advec == 'pw-scheme' ) THEN |
---|
461 | WRITE ( io, 116 ) |
---|
462 | ELSEIF ( scalar_advec == 'ws-scheme' ) THEN |
---|
463 | WRITE ( io, 504 ) |
---|
464 | ELSEIF ( scalar_advec == 'ws-scheme-mono' ) THEN |
---|
465 | WRITE ( io, 513 ) |
---|
466 | ELSE |
---|
467 | WRITE ( io, 118 ) |
---|
468 | ENDIF |
---|
469 | |
---|
470 | WRITE ( io, 139 ) TRIM( loop_optimization ) |
---|
471 | |
---|
472 | IF ( galilei_transformation ) THEN |
---|
473 | IF ( use_ug_for_galilei_tr ) THEN |
---|
474 | char1 = '0.6 * geostrophic wind' |
---|
475 | ELSE |
---|
476 | char1 = 'mean wind in model domain' |
---|
477 | ENDIF |
---|
478 | IF ( simulated_time_at_begin == simulated_time ) THEN |
---|
479 | char2 = 'at the start of the run' |
---|
480 | ELSE |
---|
481 | char2 = 'at the end of the run' |
---|
482 | ENDIF |
---|
483 | WRITE ( io, 119 ) TRIM( char1 ), TRIM( char2 ), & |
---|
484 | advected_distance_x/1000.0_wp, & |
---|
485 | advected_distance_y/1000.0_wp |
---|
486 | ENDIF |
---|
487 | WRITE ( io, 122 ) timestep_scheme |
---|
488 | IF ( use_upstream_for_tke ) WRITE ( io, 143 ) |
---|
489 | IF ( rayleigh_damping_factor /= 0.0_wp ) THEN |
---|
490 | IF ( .NOT. ocean ) THEN |
---|
491 | WRITE ( io, 123 ) 'above', rayleigh_damping_height, & |
---|
492 | rayleigh_damping_factor |
---|
493 | ELSE |
---|
494 | WRITE ( io, 123 ) 'below', rayleigh_damping_height, & |
---|
495 | rayleigh_damping_factor |
---|
496 | ENDIF |
---|
497 | ENDIF |
---|
498 | IF ( neutral ) WRITE ( io, 131 ) pt_surface |
---|
499 | IF ( humidity ) THEN |
---|
500 | IF ( .NOT. cloud_physics ) THEN |
---|
501 | WRITE ( io, 129 ) |
---|
502 | ELSE |
---|
503 | WRITE ( io, 130 ) |
---|
504 | ENDIF |
---|
505 | ENDIF |
---|
506 | IF ( passive_scalar ) WRITE ( io, 134 ) |
---|
507 | IF ( conserve_volume_flow ) THEN |
---|
508 | WRITE ( io, 150 ) conserve_volume_flow_mode |
---|
509 | IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' ) THEN |
---|
510 | WRITE ( io, 151 ) u_bulk, v_bulk |
---|
511 | ENDIF |
---|
512 | ELSEIF ( dp_external ) THEN |
---|
513 | IF ( dp_smooth ) THEN |
---|
514 | WRITE ( io, 152 ) dpdxy, dp_level_b, ', vertically smoothed.' |
---|
515 | ELSE |
---|
516 | WRITE ( io, 152 ) dpdxy, dp_level_b, '.' |
---|
517 | ENDIF |
---|
518 | ENDIF |
---|
519 | WRITE ( io, 99 ) |
---|
520 | |
---|
521 | ! |
---|
522 | !-- Runtime and timestep information |
---|
523 | WRITE ( io, 200 ) |
---|
524 | IF ( .NOT. dt_fixed ) THEN |
---|
525 | WRITE ( io, 201 ) dt_max, cfl_factor |
---|
526 | ELSE |
---|
527 | WRITE ( io, 202 ) dt |
---|
528 | ENDIF |
---|
529 | WRITE ( io, 203 ) simulated_time_at_begin, end_time |
---|
530 | |
---|
531 | IF ( time_restart /= 9999999.9_wp .AND. & |
---|
532 | simulated_time_at_begin == simulated_time ) THEN |
---|
533 | IF ( dt_restart == 9999999.9_wp ) THEN |
---|
534 | WRITE ( io, 204 ) ' Restart at: ',time_restart |
---|
535 | ELSE |
---|
536 | WRITE ( io, 205 ) ' Restart at: ',time_restart, dt_restart |
---|
537 | ENDIF |
---|
538 | ENDIF |
---|
539 | |
---|
540 | IF ( simulated_time_at_begin /= simulated_time ) THEN |
---|
541 | i = MAX ( log_point_s(10)%counts, 1 ) |
---|
542 | IF ( ( simulated_time - simulated_time_at_begin ) == 0.0_wp ) THEN |
---|
543 | cpuseconds_per_simulated_second = 0.0_wp |
---|
544 | ELSE |
---|
545 | cpuseconds_per_simulated_second = log_point_s(10)%sum / & |
---|
546 | ( simulated_time - & |
---|
547 | simulated_time_at_begin ) |
---|
548 | ENDIF |
---|
549 | WRITE ( io, 206 ) simulated_time, log_point_s(10)%sum, & |
---|
550 | log_point_s(10)%sum / REAL( i, KIND=wp ), & |
---|
551 | cpuseconds_per_simulated_second |
---|
552 | IF ( time_restart /= 9999999.9_wp .AND. time_restart < end_time ) THEN |
---|
553 | IF ( dt_restart == 9999999.9_wp ) THEN |
---|
554 | WRITE ( io, 204 ) ' Next restart at: ',time_restart |
---|
555 | ELSE |
---|
556 | WRITE ( io, 205 ) ' Next restart at: ',time_restart, dt_restart |
---|
557 | ENDIF |
---|
558 | ENDIF |
---|
559 | ENDIF |
---|
560 | |
---|
561 | |
---|
562 | ! |
---|
563 | !-- Start time for coupled runs, if independent precursor runs for atmosphere |
---|
564 | !-- and ocean are used or have been used. In this case, coupling_start_time |
---|
565 | !-- defines the time when the coupling is switched on. |
---|
566 | IF ( coupling_start_time /= 0.0_wp ) THEN |
---|
567 | WRITE ( io, 207 ) coupling_start_time |
---|
568 | ENDIF |
---|
569 | |
---|
570 | ! |
---|
571 | !-- Computational grid |
---|
572 | IF ( .NOT. ocean ) THEN |
---|
573 | WRITE ( io, 250 ) dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1) |
---|
574 | IF ( dz_stretch_level_index < nzt+1 ) THEN |
---|
575 | WRITE ( io, 252 ) dz_stretch_level, dz_stretch_level_index, & |
---|
576 | dz_stretch_factor, dz_max |
---|
577 | ENDIF |
---|
578 | ELSE |
---|
579 | WRITE ( io, 250 ) dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0) |
---|
580 | IF ( dz_stretch_level_index > 0 ) THEN |
---|
581 | WRITE ( io, 252 ) dz_stretch_level, dz_stretch_level_index, & |
---|
582 | dz_stretch_factor, dz_max |
---|
583 | ENDIF |
---|
584 | ENDIF |
---|
585 | WRITE ( io, 254 ) nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), & |
---|
586 | MIN( nnz+2, nzt+2 ) |
---|
587 | IF ( sloping_surface ) WRITE ( io, 260 ) alpha_surface |
---|
588 | |
---|
589 | ! |
---|
590 | !-- Large scale forcing and nudging |
---|
591 | WRITE ( io, 160 ) |
---|
592 | IF ( large_scale_forcing ) THEN |
---|
593 | WRITE ( io, 162 ) |
---|
594 | WRITE ( io, 163 ) |
---|
595 | |
---|
596 | IF ( large_scale_subsidence ) THEN |
---|
597 | IF ( .NOT. use_subsidence_tendencies ) THEN |
---|
598 | WRITE ( io, 164 ) |
---|
599 | ELSE |
---|
600 | WRITE ( io, 165 ) |
---|
601 | ENDIF |
---|
602 | ENDIF |
---|
603 | |
---|
604 | IF ( bc_pt_b == 'dirichlet' ) THEN |
---|
605 | WRITE ( io, 180 ) |
---|
606 | ELSEIF ( bc_pt_b == 'neumann' ) THEN |
---|
607 | WRITE ( io, 181 ) |
---|
608 | ENDIF |
---|
609 | |
---|
610 | IF ( bc_q_b == 'dirichlet' ) THEN |
---|
611 | WRITE ( io, 182 ) |
---|
612 | ELSEIF ( bc_q_b == 'neumann' ) THEN |
---|
613 | WRITE ( io, 183 ) |
---|
614 | ENDIF |
---|
615 | |
---|
616 | WRITE ( io, 167 ) |
---|
617 | IF ( nudging ) THEN |
---|
618 | WRITE ( io, 170 ) |
---|
619 | ENDIF |
---|
620 | ELSE |
---|
621 | WRITE ( io, 161 ) |
---|
622 | WRITE ( io, 171 ) |
---|
623 | ENDIF |
---|
624 | IF ( large_scale_subsidence ) THEN |
---|
625 | WRITE ( io, 168 ) |
---|
626 | WRITE ( io, 169 ) |
---|
627 | ENDIF |
---|
628 | |
---|
629 | ! |
---|
630 | !-- Profile for the large scale vertial velocity |
---|
631 | !-- Building output strings, starting with surface value |
---|
632 | IF ( large_scale_subsidence ) THEN |
---|
633 | temperatures = ' 0.0' |
---|
634 | gradients = '------' |
---|
635 | slices = ' 0' |
---|
636 | coordinates = ' 0.0' |
---|
637 | i = 1 |
---|
638 | DO WHILE ( subs_vertical_gradient_level_i(i) /= -9999 ) |
---|
639 | |
---|
640 | WRITE (coor_chr,'(E10.2,7X)') & |
---|
641 | w_subs(subs_vertical_gradient_level_i(i)) |
---|
642 | temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr ) |
---|
643 | |
---|
644 | WRITE (coor_chr,'(E10.2,7X)') subs_vertical_gradient(i) |
---|
645 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
646 | |
---|
647 | WRITE (coor_chr,'(I10,7X)') subs_vertical_gradient_level_i(i) |
---|
648 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
649 | |
---|
650 | WRITE (coor_chr,'(F10.2,7X)') subs_vertical_gradient_level(i) |
---|
651 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
652 | |
---|
653 | IF ( i == 10 ) THEN |
---|
654 | EXIT |
---|
655 | ELSE |
---|
656 | i = i + 1 |
---|
657 | ENDIF |
---|
658 | |
---|
659 | ENDDO |
---|
660 | |
---|
661 | |
---|
662 | IF ( .NOT. large_scale_forcing ) THEN |
---|
663 | WRITE ( io, 426 ) TRIM( coordinates ), TRIM( temperatures ), & |
---|
664 | TRIM( gradients ), TRIM( slices ) |
---|
665 | ENDIF |
---|
666 | |
---|
667 | |
---|
668 | ENDIF |
---|
669 | |
---|
670 | !-- Profile of the geostrophic wind (component ug) |
---|
671 | !-- Building output strings |
---|
672 | WRITE ( ugcomponent, '(F6.2)' ) ug_surface |
---|
673 | gradients = '------' |
---|
674 | slices = ' 0' |
---|
675 | coordinates = ' 0.0' |
---|
676 | i = 1 |
---|
677 | DO WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 ) |
---|
678 | |
---|
679 | WRITE (coor_chr,'(F6.2,1X)') ug(ug_vertical_gradient_level_ind(i)) |
---|
680 | ugcomponent = TRIM( ugcomponent ) // ' ' // TRIM( coor_chr ) |
---|
681 | |
---|
682 | WRITE (coor_chr,'(F6.2,1X)') ug_vertical_gradient(i) |
---|
683 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
684 | |
---|
685 | WRITE (coor_chr,'(I6,1X)') ug_vertical_gradient_level_ind(i) |
---|
686 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
687 | |
---|
688 | WRITE (coor_chr,'(F6.1,1X)') ug_vertical_gradient_level(i) |
---|
689 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
690 | |
---|
691 | IF ( i == 10 ) THEN |
---|
692 | EXIT |
---|
693 | ELSE |
---|
694 | i = i + 1 |
---|
695 | ENDIF |
---|
696 | |
---|
697 | ENDDO |
---|
698 | |
---|
699 | IF ( .NOT. large_scale_forcing ) THEN |
---|
700 | WRITE ( io, 423 ) TRIM( coordinates ), TRIM( ugcomponent ), & |
---|
701 | TRIM( gradients ), TRIM( slices ) |
---|
702 | ENDIF |
---|
703 | |
---|
704 | !-- Profile of the geostrophic wind (component vg) |
---|
705 | !-- Building output strings |
---|
706 | WRITE ( vgcomponent, '(F6.2)' ) vg_surface |
---|
707 | gradients = '------' |
---|
708 | slices = ' 0' |
---|
709 | coordinates = ' 0.0' |
---|
710 | i = 1 |
---|
711 | DO WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 ) |
---|
712 | |
---|
713 | WRITE (coor_chr,'(F6.2,1X)') vg(vg_vertical_gradient_level_ind(i)) |
---|
714 | vgcomponent = TRIM( vgcomponent ) // ' ' // TRIM( coor_chr ) |
---|
715 | |
---|
716 | WRITE (coor_chr,'(F6.2,1X)') vg_vertical_gradient(i) |
---|
717 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
718 | |
---|
719 | WRITE (coor_chr,'(I6,1X)') vg_vertical_gradient_level_ind(i) |
---|
720 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
721 | |
---|
722 | WRITE (coor_chr,'(F6.1,1X)') vg_vertical_gradient_level(i) |
---|
723 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
724 | |
---|
725 | IF ( i == 10 ) THEN |
---|
726 | EXIT |
---|
727 | ELSE |
---|
728 | i = i + 1 |
---|
729 | ENDIF |
---|
730 | |
---|
731 | ENDDO |
---|
732 | |
---|
733 | IF ( .NOT. large_scale_forcing ) THEN |
---|
734 | WRITE ( io, 424 ) TRIM( coordinates ), TRIM( vgcomponent ), & |
---|
735 | TRIM( gradients ), TRIM( slices ) |
---|
736 | ENDIF |
---|
737 | |
---|
738 | ! |
---|
739 | !-- Topography |
---|
740 | WRITE ( io, 270 ) topography |
---|
741 | SELECT CASE ( TRIM( topography ) ) |
---|
742 | |
---|
743 | CASE ( 'flat' ) |
---|
744 | ! no actions necessary |
---|
745 | |
---|
746 | CASE ( 'single_building' ) |
---|
747 | blx = INT( building_length_x / dx ) |
---|
748 | bly = INT( building_length_y / dy ) |
---|
749 | bh = INT( building_height / dz ) |
---|
750 | |
---|
751 | IF ( building_wall_left == 9999999.9_wp ) THEN |
---|
752 | building_wall_left = ( nx + 1 - blx ) / 2 * dx |
---|
753 | ENDIF |
---|
754 | bxl = INT ( building_wall_left / dx + 0.5_wp ) |
---|
755 | bxr = bxl + blx |
---|
756 | |
---|
757 | IF ( building_wall_south == 9999999.9_wp ) THEN |
---|
758 | building_wall_south = ( ny + 1 - bly ) / 2 * dy |
---|
759 | ENDIF |
---|
760 | bys = INT ( building_wall_south / dy + 0.5_wp ) |
---|
761 | byn = bys + bly |
---|
762 | |
---|
763 | WRITE ( io, 271 ) building_length_x, building_length_y, & |
---|
764 | building_height, bxl, bxr, bys, byn |
---|
765 | |
---|
766 | CASE ( 'single_street_canyon' ) |
---|
767 | ch = NINT( canyon_height / dz ) |
---|
768 | IF ( canyon_width_x /= 9999999.9_wp ) THEN |
---|
769 | ! |
---|
770 | !-- Street canyon in y direction |
---|
771 | cwx = NINT( canyon_width_x / dx ) |
---|
772 | IF ( canyon_wall_left == 9999999.9_wp ) THEN |
---|
773 | canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx |
---|
774 | ENDIF |
---|
775 | cxl = NINT( canyon_wall_left / dx ) |
---|
776 | cxr = cxl + cwx |
---|
777 | WRITE ( io, 272 ) 'y', canyon_height, ch, 'u', cxl, cxr |
---|
778 | |
---|
779 | ELSEIF ( canyon_width_y /= 9999999.9_wp ) THEN |
---|
780 | ! |
---|
781 | !-- Street canyon in x direction |
---|
782 | cwy = NINT( canyon_width_y / dy ) |
---|
783 | IF ( canyon_wall_south == 9999999.9_wp ) THEN |
---|
784 | canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy |
---|
785 | ENDIF |
---|
786 | cys = NINT( canyon_wall_south / dy ) |
---|
787 | cyn = cys + cwy |
---|
788 | WRITE ( io, 272 ) 'x', canyon_height, ch, 'v', cys, cyn |
---|
789 | ENDIF |
---|
790 | |
---|
791 | END SELECT |
---|
792 | |
---|
793 | IF ( TRIM( topography ) /= 'flat' ) THEN |
---|
794 | IF ( TRIM( topography_grid_convention ) == ' ' ) THEN |
---|
795 | IF ( TRIM( topography ) == 'single_building' .OR. & |
---|
796 | TRIM( topography ) == 'single_street_canyon' ) THEN |
---|
797 | WRITE ( io, 278 ) |
---|
798 | ELSEIF ( TRIM( topography ) == 'read_from_file' ) THEN |
---|
799 | WRITE ( io, 279 ) |
---|
800 | ENDIF |
---|
801 | ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' ) THEN |
---|
802 | WRITE ( io, 278 ) |
---|
803 | ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' ) THEN |
---|
804 | WRITE ( io, 279 ) |
---|
805 | ENDIF |
---|
806 | ENDIF |
---|
807 | |
---|
808 | IF ( plant_canopy ) THEN |
---|
809 | |
---|
810 | canopy_height = pch_index * dz |
---|
811 | |
---|
812 | WRITE ( io, 280 ) canopy_mode, canopy_height, pch_index, & |
---|
813 | canopy_drag_coeff |
---|
814 | IF ( passive_scalar ) THEN |
---|
815 | WRITE ( io, 281 ) leaf_scalar_exch_coeff, & |
---|
816 | leaf_surface_conc |
---|
817 | ENDIF |
---|
818 | |
---|
819 | ! |
---|
820 | !-- Heat flux at the top of vegetation |
---|
821 | WRITE ( io, 282 ) cthf |
---|
822 | |
---|
823 | ! |
---|
824 | !-- Leaf area density profile, calculated either from given vertical |
---|
825 | !-- gradients or from beta probability density function. |
---|
826 | IF ( .NOT. calc_beta_lad_profile ) THEN |
---|
827 | |
---|
828 | !-- Building output strings, starting with surface value |
---|
829 | WRITE ( leaf_area_density, '(F6.4)' ) lad_surface |
---|
830 | gradients = '------' |
---|
831 | slices = ' 0' |
---|
832 | coordinates = ' 0.0' |
---|
833 | i = 1 |
---|
834 | DO WHILE ( i < 11 .AND. lad_vertical_gradient_level_ind(i) /= -9999 ) |
---|
835 | |
---|
836 | WRITE (coor_chr,'(F7.2)') lad(lad_vertical_gradient_level_ind(i)) |
---|
837 | leaf_area_density = TRIM( leaf_area_density ) // ' ' // TRIM( coor_chr ) |
---|
838 | |
---|
839 | WRITE (coor_chr,'(F7.2)') lad_vertical_gradient(i) |
---|
840 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
841 | |
---|
842 | WRITE (coor_chr,'(I7)') lad_vertical_gradient_level_ind(i) |
---|
843 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
844 | |
---|
845 | WRITE (coor_chr,'(F7.1)') lad_vertical_gradient_level(i) |
---|
846 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
847 | |
---|
848 | i = i + 1 |
---|
849 | ENDDO |
---|
850 | |
---|
851 | WRITE ( io, 283 ) TRIM( coordinates ), TRIM( leaf_area_density ), & |
---|
852 | TRIM( gradients ), TRIM( slices ) |
---|
853 | |
---|
854 | ELSE |
---|
855 | |
---|
856 | WRITE ( leaf_area_density, '(F6.4)' ) lad_surface |
---|
857 | coordinates = ' 0.0' |
---|
858 | |
---|
859 | DO k = 1, pch_index |
---|
860 | |
---|
861 | WRITE (coor_chr,'(F7.2)') lad(k) |
---|
862 | leaf_area_density = TRIM( leaf_area_density ) // ' ' // TRIM( coor_chr ) |
---|
863 | |
---|
864 | WRITE (coor_chr,'(F7.1)') zu(k) |
---|
865 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
866 | |
---|
867 | ENDDO |
---|
868 | |
---|
869 | WRITE ( io, 284 ) TRIM( coordinates ), TRIM( leaf_area_density ), alpha_lad, & |
---|
870 | beta_lad, lai_beta |
---|
871 | |
---|
872 | ENDIF |
---|
873 | |
---|
874 | ENDIF |
---|
875 | |
---|
876 | |
---|
877 | IF ( land_surface ) THEN |
---|
878 | |
---|
879 | temperatures = '' |
---|
880 | gradients = '' ! use for humidity here |
---|
881 | coordinates = '' ! use for height |
---|
882 | roots = '' ! use for root fraction |
---|
883 | slices = '' ! use for index |
---|
884 | |
---|
885 | i = 1 |
---|
886 | DO i = nzb_soil, nzt_soil |
---|
887 | WRITE (coor_chr,'(F10.2,7X)') soil_temperature(i) |
---|
888 | temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr ) |
---|
889 | |
---|
890 | WRITE (coor_chr,'(F10.2,7X)') soil_moisture(i) |
---|
891 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
892 | |
---|
893 | WRITE (coor_chr,'(F10.2,7X)') - zs(i) |
---|
894 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
895 | |
---|
896 | WRITE (coor_chr,'(F10.2,7X)') root_fraction(i) |
---|
897 | roots = TRIM( roots ) // ' ' // TRIM( coor_chr ) |
---|
898 | |
---|
899 | WRITE (coor_chr,'(I10,7X)') i |
---|
900 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
901 | |
---|
902 | |
---|
903 | ENDDO |
---|
904 | |
---|
905 | ! |
---|
906 | !-- Write land surface model header |
---|
907 | WRITE( io, 419 ) |
---|
908 | IF ( conserve_water_content ) THEN |
---|
909 | WRITE( io, 440 ) |
---|
910 | ELSE |
---|
911 | WRITE( io, 441 ) |
---|
912 | ENDIF |
---|
913 | |
---|
914 | IF ( dewfall ) THEN |
---|
915 | WRITE( io, 442 ) |
---|
916 | ELSE |
---|
917 | WRITE( io, 443 ) |
---|
918 | ENDIF |
---|
919 | |
---|
920 | WRITE( io, 438 ) veg_type_name(veg_type), soil_type_name(soil_type) |
---|
921 | WRITE( io, 439 ) TRIM( coordinates ), TRIM( temperatures ), & |
---|
922 | TRIM( gradients ), TRIM( roots ), TRIM( slices ) |
---|
923 | |
---|
924 | |
---|
925 | ENDIF |
---|
926 | |
---|
927 | IF ( radiation ) THEN |
---|
928 | ! |
---|
929 | !-- Write land surface model header |
---|
930 | WRITE( io, 444 ) |
---|
931 | |
---|
932 | IF ( radiation_scheme == "constant" ) THEN |
---|
933 | WRITE( io, 445 ) net_radiation |
---|
934 | ELSEIF ( radiation_scheme == "clear-sky" ) THEN |
---|
935 | WRITE( io, 446 ) |
---|
936 | ELSE |
---|
937 | WRITE( io, 447 ) radiation_scheme |
---|
938 | ENDIF |
---|
939 | |
---|
940 | WRITE( io, 448 ) albedo |
---|
941 | WRITE( io, 449 ) dt_radiation |
---|
942 | |
---|
943 | ENDIF |
---|
944 | |
---|
945 | |
---|
946 | ! |
---|
947 | !-- Boundary conditions |
---|
948 | IF ( ibc_p_b == 0 ) THEN |
---|
949 | runten = 'p(0) = 0 |' |
---|
950 | ELSEIF ( ibc_p_b == 1 ) THEN |
---|
951 | runten = 'p(0) = p(1) |' |
---|
952 | ENDIF |
---|
953 | IF ( ibc_p_t == 0 ) THEN |
---|
954 | roben = 'p(nzt+1) = 0 |' |
---|
955 | ELSE |
---|
956 | roben = 'p(nzt+1) = p(nzt) |' |
---|
957 | ENDIF |
---|
958 | |
---|
959 | IF ( ibc_uv_b == 0 ) THEN |
---|
960 | runten = TRIM( runten ) // ' uv(0) = -uv(1) |' |
---|
961 | ELSE |
---|
962 | runten = TRIM( runten ) // ' uv(0) = uv(1) |' |
---|
963 | ENDIF |
---|
964 | IF ( TRIM( bc_uv_t ) == 'dirichlet_0' ) THEN |
---|
965 | roben = TRIM( roben ) // ' uv(nzt+1) = 0 |' |
---|
966 | ELSEIF ( ibc_uv_t == 0 ) THEN |
---|
967 | roben = TRIM( roben ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1) |' |
---|
968 | ELSE |
---|
969 | roben = TRIM( roben ) // ' uv(nzt+1) = uv(nzt) |' |
---|
970 | ENDIF |
---|
971 | |
---|
972 | IF ( ibc_pt_b == 0 ) THEN |
---|
973 | IF ( land_surface ) THEN |
---|
974 | runten = TRIM( runten ) // ' pt(0) = from soil model' |
---|
975 | ELSE |
---|
976 | runten = TRIM( runten ) // ' pt(0) = pt_surface' |
---|
977 | ENDIF |
---|
978 | ELSEIF ( ibc_pt_b == 1 ) THEN |
---|
979 | runten = TRIM( runten ) // ' pt(0) = pt(1)' |
---|
980 | ELSEIF ( ibc_pt_b == 2 ) THEN |
---|
981 | runten = TRIM( runten ) // ' pt(0) = from coupled model' |
---|
982 | ENDIF |
---|
983 | IF ( ibc_pt_t == 0 ) THEN |
---|
984 | roben = TRIM( roben ) // ' pt(nzt+1) = pt_top' |
---|
985 | ELSEIF( ibc_pt_t == 1 ) THEN |
---|
986 | roben = TRIM( roben ) // ' pt(nzt+1) = pt(nzt)' |
---|
987 | ELSEIF( ibc_pt_t == 2 ) THEN |
---|
988 | roben = TRIM( roben ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini' |
---|
989 | |
---|
990 | ENDIF |
---|
991 | |
---|
992 | WRITE ( io, 300 ) runten, roben |
---|
993 | |
---|
994 | IF ( .NOT. constant_diffusion ) THEN |
---|
995 | IF ( ibc_e_b == 1 ) THEN |
---|
996 | runten = 'e(0) = e(1)' |
---|
997 | ELSE |
---|
998 | runten = 'e(0) = e(1) = (u*/0.1)**2' |
---|
999 | ENDIF |
---|
1000 | roben = 'e(nzt+1) = e(nzt) = e(nzt-1)' |
---|
1001 | |
---|
1002 | WRITE ( io, 301 ) 'e', runten, roben |
---|
1003 | |
---|
1004 | ENDIF |
---|
1005 | |
---|
1006 | IF ( ocean ) THEN |
---|
1007 | runten = 'sa(0) = sa(1)' |
---|
1008 | IF ( ibc_sa_t == 0 ) THEN |
---|
1009 | roben = 'sa(nzt+1) = sa_surface' |
---|
1010 | ELSE |
---|
1011 | roben = 'sa(nzt+1) = sa(nzt)' |
---|
1012 | ENDIF |
---|
1013 | WRITE ( io, 301 ) 'sa', runten, roben |
---|
1014 | ENDIF |
---|
1015 | |
---|
1016 | IF ( humidity ) THEN |
---|
1017 | IF ( ibc_q_b == 0 ) THEN |
---|
1018 | IF ( land_surface ) THEN |
---|
1019 | runten = 'q(0) = from soil model' |
---|
1020 | ELSE |
---|
1021 | runten = 'q(0) = q_surface' |
---|
1022 | ENDIF |
---|
1023 | |
---|
1024 | ELSE |
---|
1025 | runten = 'q(0) = q(1)' |
---|
1026 | ENDIF |
---|
1027 | IF ( ibc_q_t == 0 ) THEN |
---|
1028 | roben = 'q(nzt) = q_top' |
---|
1029 | ELSE |
---|
1030 | roben = 'q(nzt) = q(nzt-1) + dq/dz' |
---|
1031 | ENDIF |
---|
1032 | WRITE ( io, 301 ) 'q', runten, roben |
---|
1033 | ENDIF |
---|
1034 | |
---|
1035 | IF ( passive_scalar ) THEN |
---|
1036 | IF ( ibc_q_b == 0 ) THEN |
---|
1037 | runten = 's(0) = s_surface' |
---|
1038 | ELSE |
---|
1039 | runten = 's(0) = s(1)' |
---|
1040 | ENDIF |
---|
1041 | IF ( ibc_q_t == 0 ) THEN |
---|
1042 | roben = 's(nzt) = s_top' |
---|
1043 | ELSE |
---|
1044 | roben = 's(nzt) = s(nzt-1) + ds/dz' |
---|
1045 | ENDIF |
---|
1046 | WRITE ( io, 301 ) 's', runten, roben |
---|
1047 | ENDIF |
---|
1048 | |
---|
1049 | IF ( use_surface_fluxes ) THEN |
---|
1050 | WRITE ( io, 303 ) |
---|
1051 | IF ( constant_heatflux ) THEN |
---|
1052 | IF ( large_scale_forcing .AND. lsf_surf ) THEN |
---|
1053 | WRITE ( io, 306 ) shf(0,0) |
---|
1054 | ELSE |
---|
1055 | WRITE ( io, 306 ) surface_heatflux |
---|
1056 | ENDIF |
---|
1057 | IF ( random_heatflux ) WRITE ( io, 307 ) |
---|
1058 | ENDIF |
---|
1059 | IF ( humidity .AND. constant_waterflux ) THEN |
---|
1060 | IF ( large_scale_forcing .AND. lsf_surf ) THEN |
---|
1061 | WRITE ( io, 311 ) qsws(0,0) |
---|
1062 | ELSE |
---|
1063 | WRITE ( io, 311 ) surface_waterflux |
---|
1064 | ENDIF |
---|
1065 | ENDIF |
---|
1066 | IF ( passive_scalar .AND. constant_waterflux ) THEN |
---|
1067 | WRITE ( io, 313 ) surface_waterflux |
---|
1068 | ENDIF |
---|
1069 | ENDIF |
---|
1070 | |
---|
1071 | IF ( use_top_fluxes ) THEN |
---|
1072 | WRITE ( io, 304 ) |
---|
1073 | IF ( coupling_mode == 'uncoupled' ) THEN |
---|
1074 | WRITE ( io, 320 ) top_momentumflux_u, top_momentumflux_v |
---|
1075 | IF ( constant_top_heatflux ) THEN |
---|
1076 | WRITE ( io, 306 ) top_heatflux |
---|
1077 | ENDIF |
---|
1078 | ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN |
---|
1079 | WRITE ( io, 316 ) |
---|
1080 | ENDIF |
---|
1081 | IF ( ocean .AND. constant_top_salinityflux ) THEN |
---|
1082 | WRITE ( io, 309 ) top_salinityflux |
---|
1083 | ENDIF |
---|
1084 | IF ( humidity .OR. passive_scalar ) THEN |
---|
1085 | WRITE ( io, 315 ) |
---|
1086 | ENDIF |
---|
1087 | ENDIF |
---|
1088 | |
---|
1089 | IF ( prandtl_layer ) THEN |
---|
1090 | WRITE ( io, 305 ) (zu(1)-zu(0)), roughness_length, & |
---|
1091 | z0h_factor*roughness_length, kappa, & |
---|
1092 | rif_min, rif_max |
---|
1093 | IF ( .NOT. constant_heatflux ) WRITE ( io, 308 ) |
---|
1094 | IF ( humidity .AND. .NOT. constant_waterflux ) THEN |
---|
1095 | WRITE ( io, 312 ) |
---|
1096 | ENDIF |
---|
1097 | IF ( passive_scalar .AND. .NOT. constant_waterflux ) THEN |
---|
1098 | WRITE ( io, 314 ) |
---|
1099 | ENDIF |
---|
1100 | ELSE |
---|
1101 | IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 ) THEN |
---|
1102 | WRITE ( io, 310 ) rif_min, rif_max |
---|
1103 | ENDIF |
---|
1104 | ENDIF |
---|
1105 | |
---|
1106 | WRITE ( io, 317 ) bc_lr, bc_ns |
---|
1107 | IF ( .NOT. bc_lr_cyc .OR. .NOT. bc_ns_cyc ) THEN |
---|
1108 | WRITE ( io, 318 ) use_cmax, pt_damping_width, pt_damping_factor |
---|
1109 | IF ( turbulent_inflow ) THEN |
---|
1110 | IF ( .NOT. recycling_yshift ) THEN |
---|
1111 | WRITE ( io, 319 ) recycling_width, recycling_plane, & |
---|
1112 | inflow_damping_height, inflow_damping_width |
---|
1113 | ELSE |
---|
1114 | WRITE ( io, 322 ) recycling_width, recycling_plane, & |
---|
1115 | inflow_damping_height, inflow_damping_width |
---|
1116 | END IF |
---|
1117 | ENDIF |
---|
1118 | ENDIF |
---|
1119 | |
---|
1120 | ! |
---|
1121 | !-- Initial Profiles |
---|
1122 | WRITE ( io, 321 ) |
---|
1123 | ! |
---|
1124 | !-- Initial wind profiles |
---|
1125 | IF ( u_profile(1) /= 9999999.9_wp ) WRITE ( io, 427 ) |
---|
1126 | |
---|
1127 | ! |
---|
1128 | !-- Initial temperature profile |
---|
1129 | !-- Building output strings, starting with surface temperature |
---|
1130 | WRITE ( temperatures, '(F6.2)' ) pt_surface |
---|
1131 | gradients = '------' |
---|
1132 | slices = ' 0' |
---|
1133 | coordinates = ' 0.0' |
---|
1134 | i = 1 |
---|
1135 | DO WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 ) |
---|
1136 | |
---|
1137 | WRITE (coor_chr,'(F7.2)') pt_init(pt_vertical_gradient_level_ind(i)) |
---|
1138 | temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr ) |
---|
1139 | |
---|
1140 | WRITE (coor_chr,'(F7.2)') pt_vertical_gradient(i) |
---|
1141 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
1142 | |
---|
1143 | WRITE (coor_chr,'(I7)') pt_vertical_gradient_level_ind(i) |
---|
1144 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
1145 | |
---|
1146 | WRITE (coor_chr,'(F7.1)') pt_vertical_gradient_level(i) |
---|
1147 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
1148 | |
---|
1149 | IF ( i == 10 ) THEN |
---|
1150 | EXIT |
---|
1151 | ELSE |
---|
1152 | i = i + 1 |
---|
1153 | ENDIF |
---|
1154 | |
---|
1155 | ENDDO |
---|
1156 | |
---|
1157 | IF ( .NOT. nudging ) THEN |
---|
1158 | WRITE ( io, 420 ) TRIM( coordinates ), TRIM( temperatures ), & |
---|
1159 | TRIM( gradients ), TRIM( slices ) |
---|
1160 | ELSE |
---|
1161 | WRITE ( io, 428 ) |
---|
1162 | ENDIF |
---|
1163 | |
---|
1164 | ! |
---|
1165 | !-- Initial humidity profile |
---|
1166 | !-- Building output strings, starting with surface humidity |
---|
1167 | IF ( humidity .OR. passive_scalar ) THEN |
---|
1168 | WRITE ( temperatures, '(E8.1)' ) q_surface |
---|
1169 | gradients = '--------' |
---|
1170 | slices = ' 0' |
---|
1171 | coordinates = ' 0.0' |
---|
1172 | i = 1 |
---|
1173 | DO WHILE ( q_vertical_gradient_level_ind(i) /= -9999 ) |
---|
1174 | |
---|
1175 | WRITE (coor_chr,'(E8.1,4X)') q_init(q_vertical_gradient_level_ind(i)) |
---|
1176 | temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr ) |
---|
1177 | |
---|
1178 | WRITE (coor_chr,'(E8.1,4X)') q_vertical_gradient(i) |
---|
1179 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
1180 | |
---|
1181 | WRITE (coor_chr,'(I8,4X)') q_vertical_gradient_level_ind(i) |
---|
1182 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
1183 | |
---|
1184 | WRITE (coor_chr,'(F8.1,4X)') q_vertical_gradient_level(i) |
---|
1185 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
1186 | |
---|
1187 | IF ( i == 10 ) THEN |
---|
1188 | EXIT |
---|
1189 | ELSE |
---|
1190 | i = i + 1 |
---|
1191 | ENDIF |
---|
1192 | |
---|
1193 | ENDDO |
---|
1194 | |
---|
1195 | IF ( humidity ) THEN |
---|
1196 | IF ( .NOT. nudging ) THEN |
---|
1197 | WRITE ( io, 421 ) TRIM( coordinates ), TRIM( temperatures ), & |
---|
1198 | TRIM( gradients ), TRIM( slices ) |
---|
1199 | ENDIF |
---|
1200 | ELSE |
---|
1201 | WRITE ( io, 422 ) TRIM( coordinates ), TRIM( temperatures ), & |
---|
1202 | TRIM( gradients ), TRIM( slices ) |
---|
1203 | ENDIF |
---|
1204 | ENDIF |
---|
1205 | |
---|
1206 | ! |
---|
1207 | !-- Initial salinity profile |
---|
1208 | !-- Building output strings, starting with surface salinity |
---|
1209 | IF ( ocean ) THEN |
---|
1210 | WRITE ( temperatures, '(F6.2)' ) sa_surface |
---|
1211 | gradients = '------' |
---|
1212 | slices = ' 0' |
---|
1213 | coordinates = ' 0.0' |
---|
1214 | i = 1 |
---|
1215 | DO WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 ) |
---|
1216 | |
---|
1217 | WRITE (coor_chr,'(F7.2)') sa_init(sa_vertical_gradient_level_ind(i)) |
---|
1218 | temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr ) |
---|
1219 | |
---|
1220 | WRITE (coor_chr,'(F7.2)') sa_vertical_gradient(i) |
---|
1221 | gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr ) |
---|
1222 | |
---|
1223 | WRITE (coor_chr,'(I7)') sa_vertical_gradient_level_ind(i) |
---|
1224 | slices = TRIM( slices ) // ' ' // TRIM( coor_chr ) |
---|
1225 | |
---|
1226 | WRITE (coor_chr,'(F7.1)') sa_vertical_gradient_level(i) |
---|
1227 | coordinates = TRIM( coordinates ) // ' ' // TRIM( coor_chr ) |
---|
1228 | |
---|
1229 | IF ( i == 10 ) THEN |
---|
1230 | EXIT |
---|
1231 | ELSE |
---|
1232 | i = i + 1 |
---|
1233 | ENDIF |
---|
1234 | |
---|
1235 | ENDDO |
---|
1236 | |
---|
1237 | WRITE ( io, 425 ) TRIM( coordinates ), TRIM( temperatures ), & |
---|
1238 | TRIM( gradients ), TRIM( slices ) |
---|
1239 | ENDIF |
---|
1240 | |
---|
1241 | |
---|
1242 | ! |
---|
1243 | !-- Listing of 1D-profiles |
---|
1244 | WRITE ( io, 325 ) dt_dopr_listing |
---|
1245 | IF ( averaging_interval_pr /= 0.0_wp ) THEN |
---|
1246 | WRITE ( io, 326 ) averaging_interval_pr, dt_averaging_input_pr |
---|
1247 | ENDIF |
---|
1248 | |
---|
1249 | ! |
---|
1250 | !-- DATA output |
---|
1251 | WRITE ( io, 330 ) |
---|
1252 | IF ( averaging_interval_pr /= 0.0_wp ) THEN |
---|
1253 | WRITE ( io, 326 ) averaging_interval_pr, dt_averaging_input_pr |
---|
1254 | ENDIF |
---|
1255 | |
---|
1256 | ! |
---|
1257 | !-- 1D-profiles |
---|
1258 | dopr_chr = 'Profile:' |
---|
1259 | IF ( dopr_n /= 0 ) THEN |
---|
1260 | WRITE ( io, 331 ) |
---|
1261 | |
---|
1262 | output_format = '' |
---|
1263 | output_format = output_format_netcdf |
---|
1264 | WRITE ( io, 344 ) output_format |
---|
1265 | |
---|
1266 | DO i = 1, dopr_n |
---|
1267 | dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ',' |
---|
1268 | IF ( LEN_TRIM( dopr_chr ) >= 60 ) THEN |
---|
1269 | WRITE ( io, 332 ) dopr_chr |
---|
1270 | dopr_chr = ' :' |
---|
1271 | ENDIF |
---|
1272 | ENDDO |
---|
1273 | |
---|
1274 | IF ( dopr_chr /= '' ) THEN |
---|
1275 | WRITE ( io, 332 ) dopr_chr |
---|
1276 | ENDIF |
---|
1277 | WRITE ( io, 333 ) dt_dopr, averaging_interval_pr, dt_averaging_input_pr |
---|
1278 | IF ( skip_time_dopr /= 0.0_wp ) WRITE ( io, 339 ) skip_time_dopr |
---|
1279 | ENDIF |
---|
1280 | |
---|
1281 | ! |
---|
1282 | !-- 2D-arrays |
---|
1283 | DO av = 0, 1 |
---|
1284 | |
---|
1285 | i = 1 |
---|
1286 | do2d_xy = '' |
---|
1287 | do2d_xz = '' |
---|
1288 | do2d_yz = '' |
---|
1289 | DO WHILE ( do2d(av,i) /= ' ' ) |
---|
1290 | |
---|
1291 | l = MAX( 2, LEN_TRIM( do2d(av,i) ) ) |
---|
1292 | do2d_mode = do2d(av,i)(l-1:l) |
---|
1293 | |
---|
1294 | SELECT CASE ( do2d_mode ) |
---|
1295 | CASE ( 'xy' ) |
---|
1296 | ll = LEN_TRIM( do2d_xy ) |
---|
1297 | do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ',' |
---|
1298 | CASE ( 'xz' ) |
---|
1299 | ll = LEN_TRIM( do2d_xz ) |
---|
1300 | do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ',' |
---|
1301 | CASE ( 'yz' ) |
---|
1302 | ll = LEN_TRIM( do2d_yz ) |
---|
1303 | do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ',' |
---|
1304 | END SELECT |
---|
1305 | |
---|
1306 | i = i + 1 |
---|
1307 | |
---|
1308 | ENDDO |
---|
1309 | |
---|
1310 | IF ( ( ( do2d_xy /= '' .AND. section(1,1) /= -9999 ) .OR. & |
---|
1311 | ( do2d_xz /= '' .AND. section(1,2) /= -9999 ) .OR. & |
---|
1312 | ( do2d_yz /= '' .AND. section(1,3) /= -9999 ) ) ) THEN |
---|
1313 | |
---|
1314 | IF ( av == 0 ) THEN |
---|
1315 | WRITE ( io, 334 ) '' |
---|
1316 | ELSE |
---|
1317 | WRITE ( io, 334 ) '(time-averaged)' |
---|
1318 | ENDIF |
---|
1319 | |
---|
1320 | IF ( do2d_at_begin ) THEN |
---|
1321 | begin_chr = 'and at the start' |
---|
1322 | ELSE |
---|
1323 | begin_chr = '' |
---|
1324 | ENDIF |
---|
1325 | |
---|
1326 | output_format = '' |
---|
1327 | output_format = output_format_netcdf |
---|
1328 | WRITE ( io, 344 ) output_format |
---|
1329 | |
---|
1330 | IF ( do2d_xy /= '' .AND. section(1,1) /= -9999 ) THEN |
---|
1331 | i = 1 |
---|
1332 | slices = '/' |
---|
1333 | coordinates = '/' |
---|
1334 | ! |
---|
1335 | !-- Building strings with index and coordinate information of the |
---|
1336 | !-- slices |
---|
1337 | DO WHILE ( section(i,1) /= -9999 ) |
---|
1338 | |
---|
1339 | WRITE (section_chr,'(I5)') section(i,1) |
---|
1340 | section_chr = ADJUSTL( section_chr ) |
---|
1341 | slices = TRIM( slices ) // TRIM( section_chr ) // '/' |
---|
1342 | |
---|
1343 | IF ( section(i,1) == -1 ) THEN |
---|
1344 | WRITE (coor_chr,'(F10.1)') -1.0_wp |
---|
1345 | ELSE |
---|
1346 | WRITE (coor_chr,'(F10.1)') zu(section(i,1)) |
---|
1347 | ENDIF |
---|
1348 | coor_chr = ADJUSTL( coor_chr ) |
---|
1349 | coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/' |
---|
1350 | |
---|
1351 | i = i + 1 |
---|
1352 | ENDDO |
---|
1353 | IF ( av == 0 ) THEN |
---|
1354 | WRITE ( io, 335 ) 'XY', do2d_xy, dt_do2d_xy, & |
---|
1355 | TRIM( begin_chr ), 'k', TRIM( slices ), & |
---|
1356 | TRIM( coordinates ) |
---|
1357 | IF ( skip_time_do2d_xy /= 0.0_wp ) THEN |
---|
1358 | WRITE ( io, 339 ) skip_time_do2d_xy |
---|
1359 | ENDIF |
---|
1360 | ELSE |
---|
1361 | WRITE ( io, 342 ) 'XY', do2d_xy, dt_data_output_av, & |
---|
1362 | TRIM( begin_chr ), averaging_interval, & |
---|
1363 | dt_averaging_input, 'k', TRIM( slices ), & |
---|
1364 | TRIM( coordinates ) |
---|
1365 | IF ( skip_time_data_output_av /= 0.0_wp ) THEN |
---|
1366 | WRITE ( io, 339 ) skip_time_data_output_av |
---|
1367 | ENDIF |
---|
1368 | ENDIF |
---|
1369 | IF ( netcdf_data_format > 4 ) THEN |
---|
1370 | WRITE ( io, 352 ) ntdim_2d_xy(av) |
---|
1371 | ELSE |
---|
1372 | WRITE ( io, 353 ) |
---|
1373 | ENDIF |
---|
1374 | ENDIF |
---|
1375 | |
---|
1376 | IF ( do2d_xz /= '' .AND. section(1,2) /= -9999 ) THEN |
---|
1377 | i = 1 |
---|
1378 | slices = '/' |
---|
1379 | coordinates = '/' |
---|
1380 | ! |
---|
1381 | !-- Building strings with index and coordinate information of the |
---|
1382 | !-- slices |
---|
1383 | DO WHILE ( section(i,2) /= -9999 ) |
---|
1384 | |
---|
1385 | WRITE (section_chr,'(I5)') section(i,2) |
---|
1386 | section_chr = ADJUSTL( section_chr ) |
---|
1387 | slices = TRIM( slices ) // TRIM( section_chr ) // '/' |
---|
1388 | |
---|
1389 | WRITE (coor_chr,'(F10.1)') section(i,2) * dy |
---|
1390 | coor_chr = ADJUSTL( coor_chr ) |
---|
1391 | coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/' |
---|
1392 | |
---|
1393 | i = i + 1 |
---|
1394 | ENDDO |
---|
1395 | IF ( av == 0 ) THEN |
---|
1396 | WRITE ( io, 335 ) 'XZ', do2d_xz, dt_do2d_xz, & |
---|
1397 | TRIM( begin_chr ), 'j', TRIM( slices ), & |
---|
1398 | TRIM( coordinates ) |
---|
1399 | IF ( skip_time_do2d_xz /= 0.0_wp ) THEN |
---|
1400 | WRITE ( io, 339 ) skip_time_do2d_xz |
---|
1401 | ENDIF |
---|
1402 | ELSE |
---|
1403 | WRITE ( io, 342 ) 'XZ', do2d_xz, dt_data_output_av, & |
---|
1404 | TRIM( begin_chr ), averaging_interval, & |
---|
1405 | dt_averaging_input, 'j', TRIM( slices ), & |
---|
1406 | TRIM( coordinates ) |
---|
1407 | IF ( skip_time_data_output_av /= 0.0_wp ) THEN |
---|
1408 | WRITE ( io, 339 ) skip_time_data_output_av |
---|
1409 | ENDIF |
---|
1410 | ENDIF |
---|
1411 | IF ( netcdf_data_format > 4 ) THEN |
---|
1412 | WRITE ( io, 352 ) ntdim_2d_xz(av) |
---|
1413 | ELSE |
---|
1414 | WRITE ( io, 353 ) |
---|
1415 | ENDIF |
---|
1416 | ENDIF |
---|
1417 | |
---|
1418 | IF ( do2d_yz /= '' .AND. section(1,3) /= -9999 ) THEN |
---|
1419 | i = 1 |
---|
1420 | slices = '/' |
---|
1421 | coordinates = '/' |
---|
1422 | ! |
---|
1423 | !-- Building strings with index and coordinate information of the |
---|
1424 | !-- slices |
---|
1425 | DO WHILE ( section(i,3) /= -9999 ) |
---|
1426 | |
---|
1427 | WRITE (section_chr,'(I5)') section(i,3) |
---|
1428 | section_chr = ADJUSTL( section_chr ) |
---|
1429 | slices = TRIM( slices ) // TRIM( section_chr ) // '/' |
---|
1430 | |
---|
1431 | WRITE (coor_chr,'(F10.1)') section(i,3) * dx |
---|
1432 | coor_chr = ADJUSTL( coor_chr ) |
---|
1433 | coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/' |
---|
1434 | |
---|
1435 | i = i + 1 |
---|
1436 | ENDDO |
---|
1437 | IF ( av == 0 ) THEN |
---|
1438 | WRITE ( io, 335 ) 'YZ', do2d_yz, dt_do2d_yz, & |
---|
1439 | TRIM( begin_chr ), 'i', TRIM( slices ), & |
---|
1440 | TRIM( coordinates ) |
---|
1441 | IF ( skip_time_do2d_yz /= 0.0_wp ) THEN |
---|
1442 | WRITE ( io, 339 ) skip_time_do2d_yz |
---|
1443 | ENDIF |
---|
1444 | ELSE |
---|
1445 | WRITE ( io, 342 ) 'YZ', do2d_yz, dt_data_output_av, & |
---|
1446 | TRIM( begin_chr ), averaging_interval, & |
---|
1447 | dt_averaging_input, 'i', TRIM( slices ), & |
---|
1448 | TRIM( coordinates ) |
---|
1449 | IF ( skip_time_data_output_av /= 0.0_wp ) THEN |
---|
1450 | WRITE ( io, 339 ) skip_time_data_output_av |
---|
1451 | ENDIF |
---|
1452 | ENDIF |
---|
1453 | IF ( netcdf_data_format > 4 ) THEN |
---|
1454 | WRITE ( io, 352 ) ntdim_2d_yz(av) |
---|
1455 | ELSE |
---|
1456 | WRITE ( io, 353 ) |
---|
1457 | ENDIF |
---|
1458 | ENDIF |
---|
1459 | |
---|
1460 | ENDIF |
---|
1461 | |
---|
1462 | ENDDO |
---|
1463 | |
---|
1464 | ! |
---|
1465 | !-- 3d-arrays |
---|
1466 | DO av = 0, 1 |
---|
1467 | |
---|
1468 | i = 1 |
---|
1469 | do3d_chr = '' |
---|
1470 | DO WHILE ( do3d(av,i) /= ' ' ) |
---|
1471 | |
---|
1472 | do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ',' |
---|
1473 | i = i + 1 |
---|
1474 | |
---|
1475 | ENDDO |
---|
1476 | |
---|
1477 | IF ( do3d_chr /= '' ) THEN |
---|
1478 | IF ( av == 0 ) THEN |
---|
1479 | WRITE ( io, 336 ) '' |
---|
1480 | ELSE |
---|
1481 | WRITE ( io, 336 ) '(time-averaged)' |
---|
1482 | ENDIF |
---|
1483 | |
---|
1484 | output_format = output_format_netcdf |
---|
1485 | WRITE ( io, 344 ) output_format |
---|
1486 | |
---|
1487 | IF ( do3d_at_begin ) THEN |
---|
1488 | begin_chr = 'and at the start' |
---|
1489 | ELSE |
---|
1490 | begin_chr = '' |
---|
1491 | ENDIF |
---|
1492 | IF ( av == 0 ) THEN |
---|
1493 | WRITE ( io, 337 ) do3d_chr, dt_do3d, TRIM( begin_chr ), & |
---|
1494 | zu(nz_do3d), nz_do3d |
---|
1495 | ELSE |
---|
1496 | WRITE ( io, 343 ) do3d_chr, dt_data_output_av, & |
---|
1497 | TRIM( begin_chr ), averaging_interval, & |
---|
1498 | dt_averaging_input, zu(nz_do3d), nz_do3d |
---|
1499 | ENDIF |
---|
1500 | |
---|
1501 | IF ( netcdf_data_format > 4 ) THEN |
---|
1502 | WRITE ( io, 352 ) ntdim_3d(av) |
---|
1503 | ELSE |
---|
1504 | WRITE ( io, 353 ) |
---|
1505 | ENDIF |
---|
1506 | |
---|
1507 | IF ( av == 0 ) THEN |
---|
1508 | IF ( skip_time_do3d /= 0.0_wp ) THEN |
---|
1509 | WRITE ( io, 339 ) skip_time_do3d |
---|
1510 | ENDIF |
---|
1511 | ELSE |
---|
1512 | IF ( skip_time_data_output_av /= 0.0_wp ) THEN |
---|
1513 | WRITE ( io, 339 ) skip_time_data_output_av |
---|
1514 | ENDIF |
---|
1515 | ENDIF |
---|
1516 | |
---|
1517 | ENDIF |
---|
1518 | |
---|
1519 | ENDDO |
---|
1520 | |
---|
1521 | ! |
---|
1522 | !-- masked arrays |
---|
1523 | IF ( masks > 0 ) WRITE ( io, 345 ) & |
---|
1524 | mask_scale_x, mask_scale_y, mask_scale_z |
---|
1525 | DO mid = 1, masks |
---|
1526 | DO av = 0, 1 |
---|
1527 | |
---|
1528 | i = 1 |
---|
1529 | domask_chr = '' |
---|
1530 | DO WHILE ( domask(mid,av,i) /= ' ' ) |
---|
1531 | domask_chr = TRIM( domask_chr ) // ' ' // & |
---|
1532 | TRIM( domask(mid,av,i) ) // ',' |
---|
1533 | i = i + 1 |
---|
1534 | ENDDO |
---|
1535 | |
---|
1536 | IF ( domask_chr /= '' ) THEN |
---|
1537 | IF ( av == 0 ) THEN |
---|
1538 | WRITE ( io, 346 ) '', mid |
---|
1539 | ELSE |
---|
1540 | WRITE ( io, 346 ) ' (time-averaged)', mid |
---|
1541 | ENDIF |
---|
1542 | |
---|
1543 | output_format = output_format_netcdf |
---|
1544 | !-- Parallel output not implemented for mask data, hence |
---|
1545 | !-- output_format must be adjusted. |
---|
1546 | IF ( netcdf_data_format == 5 ) output_format = 'netCDF4/HDF5' |
---|
1547 | IF ( netcdf_data_format == 6 ) output_format = 'netCDF4/HDF5 classic' |
---|
1548 | WRITE ( io, 344 ) output_format |
---|
1549 | |
---|
1550 | IF ( av == 0 ) THEN |
---|
1551 | WRITE ( io, 347 ) domask_chr, dt_domask(mid) |
---|
1552 | ELSE |
---|
1553 | WRITE ( io, 348 ) domask_chr, dt_data_output_av, & |
---|
1554 | averaging_interval, dt_averaging_input |
---|
1555 | ENDIF |
---|
1556 | |
---|
1557 | IF ( av == 0 ) THEN |
---|
1558 | IF ( skip_time_domask(mid) /= 0.0_wp ) THEN |
---|
1559 | WRITE ( io, 339 ) skip_time_domask(mid) |
---|
1560 | ENDIF |
---|
1561 | ELSE |
---|
1562 | IF ( skip_time_data_output_av /= 0.0_wp ) THEN |
---|
1563 | WRITE ( io, 339 ) skip_time_data_output_av |
---|
1564 | ENDIF |
---|
1565 | ENDIF |
---|
1566 | ! |
---|
1567 | !-- output locations |
---|
1568 | DO dim = 1, 3 |
---|
1569 | IF ( mask(mid,dim,1) >= 0.0_wp ) THEN |
---|
1570 | count = 0 |
---|
1571 | DO WHILE ( mask(mid,dim,count+1) >= 0.0_wp ) |
---|
1572 | count = count + 1 |
---|
1573 | ENDDO |
---|
1574 | WRITE ( io, 349 ) dir(dim), dir(dim), mid, dir(dim), & |
---|
1575 | mask(mid,dim,:count) |
---|
1576 | ELSEIF ( mask_loop(mid,dim,1) < 0.0_wp .AND. & |
---|
1577 | mask_loop(mid,dim,2) < 0.0_wp .AND. & |
---|
1578 | mask_loop(mid,dim,3) == 0.0_wp ) THEN |
---|
1579 | WRITE ( io, 350 ) dir(dim), dir(dim) |
---|
1580 | ELSEIF ( mask_loop(mid,dim,3) == 0.0_wp ) THEN |
---|
1581 | WRITE ( io, 351 ) dir(dim), dir(dim), mid, dir(dim), & |
---|
1582 | mask_loop(mid,dim,1:2) |
---|
1583 | ELSE |
---|
1584 | WRITE ( io, 351 ) dir(dim), dir(dim), mid, dir(dim), & |
---|
1585 | mask_loop(mid,dim,1:3) |
---|
1586 | ENDIF |
---|
1587 | ENDDO |
---|
1588 | ENDIF |
---|
1589 | |
---|
1590 | ENDDO |
---|
1591 | ENDDO |
---|
1592 | |
---|
1593 | ! |
---|
1594 | !-- Timeseries |
---|
1595 | IF ( dt_dots /= 9999999.9_wp ) THEN |
---|
1596 | WRITE ( io, 340 ) |
---|
1597 | |
---|
1598 | output_format = output_format_netcdf |
---|
1599 | WRITE ( io, 344 ) output_format |
---|
1600 | WRITE ( io, 341 ) dt_dots |
---|
1601 | ENDIF |
---|
1602 | |
---|
1603 | #if defined( __dvrp_graphics ) |
---|
1604 | ! |
---|
1605 | !-- Dvrp-output |
---|
1606 | IF ( dt_dvrp /= 9999999.9_wp ) THEN |
---|
1607 | WRITE ( io, 360 ) dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), & |
---|
1608 | TRIM( dvrp_username ), TRIM( dvrp_directory ) |
---|
1609 | i = 1 |
---|
1610 | l = 0 |
---|
1611 | m = 0 |
---|
1612 | DO WHILE ( mode_dvrp(i) /= ' ' ) |
---|
1613 | IF ( mode_dvrp(i)(1:10) == 'isosurface' ) THEN |
---|
1614 | READ ( mode_dvrp(i), '(10X,I2)' ) j |
---|
1615 | l = l + 1 |
---|
1616 | IF ( do3d(0,j) /= ' ' ) THEN |
---|
1617 | WRITE ( io, 361 ) TRIM( do3d(0,j) ), threshold(l), & |
---|
1618 | isosurface_color(:,l) |
---|
1619 | ENDIF |
---|
1620 | ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' ) THEN |
---|
1621 | READ ( mode_dvrp(i), '(6X,I2)' ) j |
---|
1622 | m = m + 1 |
---|
1623 | IF ( do2d(0,j) /= ' ' ) THEN |
---|
1624 | WRITE ( io, 362 ) TRIM( do2d(0,j) ), & |
---|
1625 | slicer_range_limits_dvrp(:,m) |
---|
1626 | ENDIF |
---|
1627 | ELSEIF ( mode_dvrp(i)(1:9) == 'particles' ) THEN |
---|
1628 | WRITE ( io, 363 ) dvrp_psize |
---|
1629 | IF ( particle_dvrpsize /= 'none' ) THEN |
---|
1630 | WRITE ( io, 364 ) 'size', TRIM( particle_dvrpsize ), & |
---|
1631 | dvrpsize_interval |
---|
1632 | ENDIF |
---|
1633 | IF ( particle_color /= 'none' ) THEN |
---|
1634 | WRITE ( io, 364 ) 'color', TRIM( particle_color ), & |
---|
1635 | color_interval |
---|
1636 | ENDIF |
---|
1637 | ENDIF |
---|
1638 | i = i + 1 |
---|
1639 | ENDDO |
---|
1640 | |
---|
1641 | WRITE ( io, 365 ) groundplate_color, superelevation_x, & |
---|
1642 | superelevation_y, superelevation, clip_dvrp_l, & |
---|
1643 | clip_dvrp_r, clip_dvrp_s, clip_dvrp_n |
---|
1644 | |
---|
1645 | IF ( TRIM( topography ) /= 'flat' ) THEN |
---|
1646 | WRITE ( io, 366 ) topography_color |
---|
1647 | IF ( cluster_size > 1 ) THEN |
---|
1648 | WRITE ( io, 367 ) cluster_size |
---|
1649 | ENDIF |
---|
1650 | ENDIF |
---|
1651 | |
---|
1652 | ENDIF |
---|
1653 | #endif |
---|
1654 | |
---|
1655 | #if defined( __spectra ) |
---|
1656 | ! |
---|
1657 | !-- Spectra output |
---|
1658 | IF ( dt_dosp /= 9999999.9_wp ) THEN |
---|
1659 | WRITE ( io, 370 ) |
---|
1660 | |
---|
1661 | output_format = output_format_netcdf |
---|
1662 | WRITE ( io, 344 ) output_format |
---|
1663 | WRITE ( io, 371 ) dt_dosp |
---|
1664 | IF ( skip_time_dosp /= 0.0_wp ) WRITE ( io, 339 ) skip_time_dosp |
---|
1665 | WRITE ( io, 372 ) ( data_output_sp(i), i = 1,10 ), & |
---|
1666 | ( spectra_direction(i), i = 1,10 ), & |
---|
1667 | ( comp_spectra_level(i), i = 1,100 ), & |
---|
1668 | ( plot_spectra_level(i), i = 1,100 ), & |
---|
1669 | averaging_interval_sp, dt_averaging_input_pr |
---|
1670 | ENDIF |
---|
1671 | #endif |
---|
1672 | |
---|
1673 | WRITE ( io, 99 ) |
---|
1674 | |
---|
1675 | ! |
---|
1676 | !-- Physical quantities |
---|
1677 | WRITE ( io, 400 ) |
---|
1678 | |
---|
1679 | ! |
---|
1680 | !-- Geostrophic parameters |
---|
1681 | IF ( radiation .AND. radiation_scheme /= 'constant' ) THEN |
---|
1682 | WRITE ( io, 417 ) lambda |
---|
1683 | ENDIF |
---|
1684 | WRITE ( io, 410 ) phi, omega, f, fs |
---|
1685 | |
---|
1686 | ! |
---|
1687 | !-- Other quantities |
---|
1688 | WRITE ( io, 411 ) g |
---|
1689 | IF ( radiation .AND. radiation_scheme /= 'constant' ) THEN |
---|
1690 | WRITE ( io, 418 ) day_init, time_utc_init |
---|
1691 | ENDIF |
---|
1692 | |
---|
1693 | WRITE ( io, 412 ) TRIM( reference_state ) |
---|
1694 | IF ( use_single_reference_value ) THEN |
---|
1695 | IF ( ocean ) THEN |
---|
1696 | WRITE ( io, 413 ) prho_reference |
---|
1697 | ELSE |
---|
1698 | WRITE ( io, 414 ) pt_reference |
---|
1699 | ENDIF |
---|
1700 | ENDIF |
---|
1701 | |
---|
1702 | ! |
---|
1703 | !-- Cloud physics parameters |
---|
1704 | IF ( cloud_physics ) THEN |
---|
1705 | WRITE ( io, 415 ) |
---|
1706 | WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v |
---|
1707 | IF ( icloud_scheme == 0 ) THEN |
---|
1708 | WRITE ( io, 510 ) 1.0E-6_wp * nc_const |
---|
1709 | IF ( precipitation ) WRITE ( io, 511 ) c_sedimentation |
---|
1710 | ENDIF |
---|
1711 | ENDIF |
---|
1712 | |
---|
1713 | ! |
---|
1714 | !-- Cloud physcis parameters / quantities / numerical methods |
---|
1715 | WRITE ( io, 430 ) |
---|
1716 | IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets) THEN |
---|
1717 | WRITE ( io, 431 ) |
---|
1718 | ELSEIF ( humidity .AND. cloud_physics ) THEN |
---|
1719 | WRITE ( io, 432 ) |
---|
1720 | IF ( cloud_top_radiation ) WRITE ( io, 132 ) |
---|
1721 | IF ( icloud_scheme == 1 ) THEN |
---|
1722 | IF ( precipitation ) WRITE ( io, 133 ) |
---|
1723 | ELSEIF ( icloud_scheme == 0 ) THEN |
---|
1724 | IF ( drizzle ) WRITE ( io, 506 ) |
---|
1725 | IF ( precipitation ) THEN |
---|
1726 | WRITE ( io, 505 ) |
---|
1727 | IF ( turbulence ) WRITE ( io, 507 ) |
---|
1728 | IF ( ventilation_effect ) WRITE ( io, 508 ) |
---|
1729 | IF ( limiter_sedimentation ) WRITE ( io, 509 ) |
---|
1730 | ENDIF |
---|
1731 | ENDIF |
---|
1732 | ELSEIF ( humidity .AND. cloud_droplets ) THEN |
---|
1733 | WRITE ( io, 433 ) |
---|
1734 | IF ( curvature_solution_effects ) WRITE ( io, 434 ) |
---|
1735 | IF ( collision_kernel /= 'none' ) THEN |
---|
1736 | WRITE ( io, 435 ) TRIM( collision_kernel ) |
---|
1737 | IF ( collision_kernel(6:9) == 'fast' ) THEN |
---|
1738 | WRITE ( io, 436 ) radius_classes, dissipation_classes |
---|
1739 | ENDIF |
---|
1740 | ELSE |
---|
1741 | WRITE ( io, 437 ) |
---|
1742 | ENDIF |
---|
1743 | ENDIF |
---|
1744 | |
---|
1745 | ! |
---|
1746 | !-- LES / turbulence parameters |
---|
1747 | WRITE ( io, 450 ) |
---|
1748 | |
---|
1749 | !-- |
---|
1750 | ! ... LES-constants used must still be added here |
---|
1751 | !-- |
---|
1752 | IF ( constant_diffusion ) THEN |
---|
1753 | WRITE ( io, 451 ) km_constant, km_constant/prandtl_number, & |
---|
1754 | prandtl_number |
---|
1755 | ENDIF |
---|
1756 | IF ( .NOT. constant_diffusion) THEN |
---|
1757 | IF ( e_init > 0.0_wp ) WRITE ( io, 455 ) e_init |
---|
1758 | IF ( e_min > 0.0_wp ) WRITE ( io, 454 ) e_min |
---|
1759 | IF ( wall_adjustment ) WRITE ( io, 453 ) wall_adjustment_factor |
---|
1760 | ENDIF |
---|
1761 | |
---|
1762 | ! |
---|
1763 | !-- Special actions during the run |
---|
1764 | WRITE ( io, 470 ) |
---|
1765 | IF ( create_disturbances ) THEN |
---|
1766 | WRITE ( io, 471 ) dt_disturb, disturbance_amplitude, & |
---|
1767 | zu(disturbance_level_ind_b), disturbance_level_ind_b,& |
---|
1768 | zu(disturbance_level_ind_t), disturbance_level_ind_t |
---|
1769 | IF ( .NOT. bc_lr_cyc .OR. .NOT. bc_ns_cyc ) THEN |
---|
1770 | WRITE ( io, 472 ) inflow_disturbance_begin, inflow_disturbance_end |
---|
1771 | ELSE |
---|
1772 | WRITE ( io, 473 ) disturbance_energy_limit |
---|
1773 | ENDIF |
---|
1774 | WRITE ( io, 474 ) TRIM( random_generator ) |
---|
1775 | ENDIF |
---|
1776 | IF ( pt_surface_initial_change /= 0.0_wp ) THEN |
---|
1777 | WRITE ( io, 475 ) pt_surface_initial_change |
---|
1778 | ENDIF |
---|
1779 | IF ( humidity .AND. q_surface_initial_change /= 0.0_wp ) THEN |
---|
1780 | WRITE ( io, 476 ) q_surface_initial_change |
---|
1781 | ENDIF |
---|
1782 | IF ( passive_scalar .AND. q_surface_initial_change /= 0.0_wp ) THEN |
---|
1783 | WRITE ( io, 477 ) q_surface_initial_change |
---|
1784 | ENDIF |
---|
1785 | |
---|
1786 | IF ( particle_advection ) THEN |
---|
1787 | ! |
---|
1788 | !-- Particle attributes |
---|
1789 | WRITE ( io, 480 ) particle_advection_start, dt_prel, bc_par_lr, & |
---|
1790 | bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, & |
---|
1791 | end_time_prel |
---|
1792 | IF ( use_sgs_for_particles ) WRITE ( io, 488 ) dt_min_part |
---|
1793 | IF ( random_start_position ) WRITE ( io, 481 ) |
---|
1794 | IF ( seed_follows_topography ) WRITE ( io, 496 ) |
---|
1795 | IF ( particles_per_point > 1 ) WRITE ( io, 489 ) particles_per_point |
---|
1796 | WRITE ( io, 495 ) total_number_of_particles |
---|
1797 | IF ( use_particle_tails .AND. maximum_number_of_tailpoints /= 0 ) THEN |
---|
1798 | WRITE ( io, 483 ) maximum_number_of_tailpoints |
---|
1799 | IF ( minimum_tailpoint_distance /= 0 ) THEN |
---|
1800 | WRITE ( io, 484 ) total_number_of_tails, & |
---|
1801 | minimum_tailpoint_distance, & |
---|
1802 | maximum_tailpoint_age |
---|
1803 | ENDIF |
---|
1804 | ENDIF |
---|
1805 | IF ( dt_write_particle_data /= 9999999.9_wp ) THEN |
---|
1806 | WRITE ( io, 485 ) dt_write_particle_data |
---|
1807 | IF ( netcdf_data_format > 1 ) THEN |
---|
1808 | output_format = 'netcdf (64 bit offset) and binary' |
---|
1809 | ELSE |
---|
1810 | output_format = 'netcdf and binary' |
---|
1811 | ENDIF |
---|
1812 | WRITE ( io, 344 ) output_format |
---|
1813 | ENDIF |
---|
1814 | IF ( dt_dopts /= 9999999.9_wp ) WRITE ( io, 494 ) dt_dopts |
---|
1815 | IF ( write_particle_statistics ) WRITE ( io, 486 ) |
---|
1816 | |
---|
1817 | WRITE ( io, 487 ) number_of_particle_groups |
---|
1818 | |
---|
1819 | DO i = 1, number_of_particle_groups |
---|
1820 | IF ( i == 1 .AND. density_ratio(i) == 9999999.9_wp ) THEN |
---|
1821 | WRITE ( io, 490 ) i, 0.0_wp |
---|
1822 | WRITE ( io, 492 ) |
---|
1823 | ELSE |
---|
1824 | WRITE ( io, 490 ) i, radius(i) |
---|
1825 | IF ( density_ratio(i) /= 0.0_wp ) THEN |
---|
1826 | WRITE ( io, 491 ) density_ratio(i) |
---|
1827 | ELSE |
---|
1828 | WRITE ( io, 492 ) |
---|
1829 | ENDIF |
---|
1830 | ENDIF |
---|
1831 | WRITE ( io, 493 ) psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), & |
---|
1832 | pdx(i), pdy(i), pdz(i) |
---|
1833 | IF ( .NOT. vertical_particle_advection(i) ) WRITE ( io, 482 ) |
---|
1834 | ENDDO |
---|
1835 | |
---|
1836 | ENDIF |
---|
1837 | |
---|
1838 | |
---|
1839 | ! |
---|
1840 | !-- Parameters of 1D-model |
---|
1841 | IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 ) THEN |
---|
1842 | WRITE ( io, 500 ) end_time_1d, dt_run_control_1d, dt_pr_1d, & |
---|
1843 | mixing_length_1d, dissipation_1d |
---|
1844 | IF ( damp_level_ind_1d /= nzt+1 ) THEN |
---|
1845 | WRITE ( io, 502 ) zu(damp_level_ind_1d), damp_level_ind_1d |
---|
1846 | ENDIF |
---|
1847 | ENDIF |
---|
1848 | |
---|
1849 | ! |
---|
1850 | !-- User-defined information |
---|
1851 | CALL user_header( io ) |
---|
1852 | |
---|
1853 | WRITE ( io, 99 ) |
---|
1854 | |
---|
1855 | ! |
---|
1856 | !-- Write buffer contents to disc immediately |
---|
1857 | CALL local_flush( io ) |
---|
1858 | |
---|
1859 | ! |
---|
1860 | !-- Here the FORMATs start |
---|
1861 | |
---|
1862 | 99 FORMAT (1X,78('-')) |
---|
1863 | 100 FORMAT (/1X,'******************************',4X,44('-')/ & |
---|
1864 | 1X,'* ',A,' *',4X,A/ & |
---|
1865 | 1X,'******************************',4X,44('-')) |
---|
1866 | 101 FORMAT (35X,'coupled run using MPI-',I1,': ',A/ & |
---|
1867 | 35X,42('-')) |
---|
1868 | 102 FORMAT (/' Date: ',A8,4X,'Run: ',A20/ & |
---|
1869 | ' Time: ',A8,4X,'Run-No.: ',I2.2/ & |
---|
1870 | ' Run on host: ',A10) |
---|
1871 | #if defined( __parallel ) |
---|
1872 | 103 FORMAT (' Number of PEs:',10X,I6,4X,'Processor grid (x,y): (',I4,',',I4, & |
---|
1873 | ')',1X,A) |
---|
1874 | 104 FORMAT (' Number of PEs:',10X,I6,4X,'Tasks:',I4,' threads per task:',I4/ & |
---|
1875 | 35X,'Processor grid (x,y): (',I4,',',I4,')',1X,A) |
---|
1876 | 105 FORMAT (35X,'One additional PE is used to handle'/37X,'the dvrp output!') |
---|
1877 | 106 FORMAT (35X,'A 1d-decomposition along x is forced'/ & |
---|
1878 | 35X,'because the job is running on an SMP-cluster') |
---|
1879 | 107 FORMAT (35X,'A 1d-decomposition along ',A,' is used') |
---|
1880 | 108 FORMAT (35X,'Max. # of parallel I/O streams is ',I5) |
---|
1881 | 109 FORMAT (35X,'Precursor run for coupled atmos-ocean run'/ & |
---|
1882 | 35X,42('-')) |
---|
1883 | 114 FORMAT (35X,'Coupled atmosphere-ocean run following'/ & |
---|
1884 | 35X,'independent precursor runs'/ & |
---|
1885 | 35X,42('-')) |
---|
1886 | 117 FORMAT (' Accelerator boards / node: ',I2) |
---|
1887 | #endif |
---|
1888 | 110 FORMAT (/' Numerical Schemes:'/ & |
---|
1889 | ' -----------------'/) |
---|
1890 | 111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines') |
---|
1891 | 112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ & |
---|
1892 | ' Iterations (initial/other): ',I3,'/',I3,' omega = ',F5.3) |
---|
1893 | 113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', & |
---|
1894 | ' or Upstream') |
---|
1895 | 115 FORMAT (' FFT and transpositions are overlapping') |
---|
1896 | 116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', & |
---|
1897 | ' or Upstream') |
---|
1898 | 118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme') |
---|
1899 | 119 FORMAT (' --> Galilei-Transform applied to horizontal advection:'/ & |
---|
1900 | ' translation velocity = ',A/ & |
---|
1901 | ' distance advected ',A,': ',F8.3,' km(x) ',F8.3,' km(y)') |
---|
1902 | 120 FORMAT (' Accelerator boards: ',8X,I2) |
---|
1903 | 122 FORMAT (' --> Time differencing scheme: ',A) |
---|
1904 | 123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ & |
---|
1905 | ' maximum damping coefficient: ',F5.3, ' 1/s') |
---|
1906 | 129 FORMAT (' --> Additional prognostic equation for the specific humidity') |
---|
1907 | 130 FORMAT (' --> Additional prognostic equation for the total water content') |
---|
1908 | 131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', & |
---|
1909 | F6.2, ' K assumed') |
---|
1910 | 132 FORMAT (' Parameterization of long-wave radiation processes via'/ & |
---|
1911 | ' effective emissivity scheme') |
---|
1912 | 133 FORMAT (' Precipitation parameterization via Kessler-Scheme') |
---|
1913 | 134 FORMAT (' --> Additional prognostic equation for a passive scalar') |
---|
1914 | 135 FORMAT (' --> Solve perturbation pressure via ',A,' method (', & |
---|
1915 | A,'-cycle)'/ & |
---|
1916 | ' number of grid levels: ',I2/ & |
---|
1917 | ' Gauss-Seidel red/black iterations: ',I2) |
---|
1918 | 136 FORMAT (' gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', & |
---|
1919 | I3,')') |
---|
1920 | 137 FORMAT (' level data gathered on PE0 at level: ',I2/ & |
---|
1921 | ' gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', & |
---|
1922 | I3,')'/ & |
---|
1923 | ' gridpoints of coarsest domain (x,y,z): (',I3,',',I3,',', & |
---|
1924 | I3,')') |
---|
1925 | 139 FORMAT (' --> Loop optimization method: ',A) |
---|
1926 | 140 FORMAT (' maximum residual allowed: ',E10.3) |
---|
1927 | 141 FORMAT (' fixed number of multigrid cycles: ',I4) |
---|
1928 | 142 FORMAT (' perturbation pressure is calculated at every Runge-Kutta ', & |
---|
1929 | 'step') |
---|
1930 | 143 FORMAT (' Euler/upstream scheme is used for the SGS turbulent ', & |
---|
1931 | 'kinetic energy') |
---|
1932 | 144 FORMAT (' masking method is used') |
---|
1933 | 150 FORMAT (' --> Volume flow at the right and north boundary will be ', & |
---|
1934 | 'conserved'/ & |
---|
1935 | ' using the ',A,' mode') |
---|
1936 | 151 FORMAT (' with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s') |
---|
1937 | 152 FORMAT (' --> External pressure gradient directly prescribed by the user:',& |
---|
1938 | /' ',2(1X,E12.5),'Pa/m in x/y direction', & |
---|
1939 | /' starting from dp_level_b =', F8.3, 'm', A /) |
---|
1940 | 160 FORMAT (//' Large scale forcing and nudging:'/ & |
---|
1941 | ' -------------------------------'/) |
---|
1942 | 161 FORMAT (' --> No large scale forcing from external is used (default) ') |
---|
1943 | 162 FORMAT (' --> Large scale forcing from external file LSF_DATA is used: ') |
---|
1944 | 163 FORMAT (' - large scale advection tendencies ') |
---|
1945 | 164 FORMAT (' - large scale subsidence velocity w_subs ') |
---|
1946 | 165 FORMAT (' - large scale subsidence tendencies ') |
---|
1947 | 167 FORMAT (' - and geostrophic wind components ug and vg') |
---|
1948 | 168 FORMAT (' --> Large-scale vertical motion is used in the ', & |
---|
1949 | 'prognostic equation(s) for') |
---|
1950 | 169 FORMAT (' the scalar(s) only') |
---|
1951 | 170 FORMAT (' --> Nudging is used') |
---|
1952 | 171 FORMAT (' --> No nudging is used (default) ') |
---|
1953 | 180 FORMAT (' - prescribed surface values for temperature') |
---|
1954 | 181 FORMAT (' - prescribed surface fluxes for temperature') |
---|
1955 | 182 FORMAT (' - prescribed surface values for humidity') |
---|
1956 | 183 FORMAT (' - prescribed surface fluxes for humidity') |
---|
1957 | 200 FORMAT (//' Run time and time step information:'/ & |
---|
1958 | ' ----------------------------------'/) |
---|
1959 | 201 FORMAT ( ' Timestep: variable maximum value: ',F6.3,' s', & |
---|
1960 | ' CFL-factor: ',F4.2) |
---|
1961 | 202 FORMAT ( ' Timestep: dt = ',F6.3,' s'/) |
---|
1962 | 203 FORMAT ( ' Start time: ',F9.3,' s'/ & |
---|
1963 | ' End time: ',F9.3,' s') |
---|
1964 | 204 FORMAT ( A,F9.3,' s') |
---|
1965 | 205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s') |
---|
1966 | 206 FORMAT (/' Time reached: ',F9.3,' s'/ & |
---|
1967 | ' CPU-time used: ',F9.3,' s per timestep: ', & |
---|
1968 | ' ',F9.3,' s'/ & |
---|
1969 | ' per second of simulated tim', & |
---|
1970 | 'e: ',F9.3,' s') |
---|
1971 | 207 FORMAT ( ' Coupling start time: ',F9.3,' s') |
---|
1972 | 250 FORMAT (//' Computational grid and domain size:'/ & |
---|
1973 | ' ----------------------------------'// & |
---|
1974 | ' Grid length: dx = ',F7.3,' m dy = ',F7.3, & |
---|
1975 | ' m dz = ',F7.3,' m'/ & |
---|
1976 | ' Domain size: x = ',F10.3,' m y = ',F10.3, & |
---|
1977 | ' m z(u) = ',F10.3,' m'/) |
---|
1978 | 252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', & |
---|
1979 | ' factor: ',F5.3/ & |
---|
1980 | ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/) |
---|
1981 | 254 FORMAT (' Number of gridpoints (x,y,z): (0:',I4,', 0:',I4,', 0:',I4,')'/ & |
---|
1982 | ' Subdomain size (x,y,z): ( ',I4,', ',I4,', ',I4,')'/) |
---|
1983 | 260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,& |
---|
1984 | ' degrees') |
---|
1985 | 270 FORMAT (//' Topography information:'/ & |
---|
1986 | ' ----------------------'// & |
---|
1987 | 1X,'Topography: ',A) |
---|
1988 | 271 FORMAT ( ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ & |
---|
1989 | ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, & |
---|
1990 | ' / ',I4) |
---|
1991 | 272 FORMAT ( ' Single quasi-2D street canyon of infinite length in ',A, & |
---|
1992 | ' direction' / & |
---|
1993 | ' Canyon height: ', F6.2, 'm, ch = ', I4, '.' / & |
---|
1994 | ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.') |
---|
1995 | 278 FORMAT (' Topography grid definition convention:'/ & |
---|
1996 | ' cell edge (staggered grid points'/ & |
---|
1997 | ' (u in x-direction, v in y-direction))' /) |
---|
1998 | 279 FORMAT (' Topography grid definition convention:'/ & |
---|
1999 | ' cell center (scalar grid points)' /) |
---|
2000 | 280 FORMAT (//' Vegetation canopy (drag) model:'/ & |
---|
2001 | ' ------------------------------'// & |
---|
2002 | ' Canopy mode: ', A / & |
---|
2003 | ' Canopy height: ',F6.2,'m (',I4,' grid points)' / & |
---|
2004 | ' Leaf drag coefficient: ',F6.2 /) |
---|
2005 | 281 FORMAT (/ ' Scalar exchange coefficient: ',F6.2 / & |
---|
2006 | ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /) |
---|
2007 | 282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s') |
---|
2008 | 283 FORMAT (/ ' Characteristic levels of the leaf area density:'// & |
---|
2009 | ' Height: ',A,' m'/ & |
---|
2010 | ' Leaf area density: ',A,' m**2/m**3'/ & |
---|
2011 | ' Gradient: ',A,' m**2/m**4'/ & |
---|
2012 | ' Gridpoint: ',A) |
---|
2013 | 284 FORMAT (//' Characteristic levels of the leaf area density and coefficients:'// & |
---|
2014 | ' Height: ',A,' m'/ & |
---|
2015 | ' Leaf area density: ',A,' m**2/m**3'/ & |
---|
2016 | ' Coefficient alpha: ',F6.2 / & |
---|
2017 | ' Coefficient beta: ',F6.2 / & |
---|
2018 | ' Leaf area index: ',F6.2,' m**2/m**2' /) |
---|
2019 | |
---|
2020 | 300 FORMAT (//' Boundary conditions:'/ & |
---|
2021 | ' -------------------'// & |
---|
2022 | ' p uv ', & |
---|
2023 | ' pt'// & |
---|
2024 | ' B. bound.: ',A/ & |
---|
2025 | ' T. bound.: ',A) |
---|
2026 | 301 FORMAT (/' ',A// & |
---|
2027 | ' B. bound.: ',A/ & |
---|
2028 | ' T. bound.: ',A) |
---|
2029 | 303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1') |
---|
2030 | 304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt') |
---|
2031 | 305 FORMAT (//' Prandtl-Layer between bottom surface and first ', & |
---|
2032 | 'computational u,v-level:'// & |
---|
2033 | ' zp = ',F6.2,' m z0 = ',F6.4,' m z0h = ',F7.5,& |
---|
2034 | ' m kappa = ',F4.2/ & |
---|
2035 | ' Rif value range: ',F6.2,' <= rif <=',F6.2) |
---|
2036 | 306 FORMAT (' Predefined constant heatflux: ',F9.6,' K m/s') |
---|
2037 | 307 FORMAT (' Heatflux has a random normal distribution') |
---|
2038 | 308 FORMAT (' Predefined surface temperature') |
---|
2039 | 309 FORMAT (' Predefined constant salinityflux: ',F9.6,' psu m/s') |
---|
2040 | 310 FORMAT (//' 1D-Model:'// & |
---|
2041 | ' Rif value range: ',F6.2,' <= rif <=',F6.2) |
---|
2042 | 311 FORMAT (' Predefined constant humidity flux: ',E10.3,' m/s') |
---|
2043 | 312 FORMAT (' Predefined surface humidity') |
---|
2044 | 313 FORMAT (' Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)') |
---|
2045 | 314 FORMAT (' Predefined scalar value at the surface') |
---|
2046 | 315 FORMAT (' Humidity / scalar flux at top surface is 0.0') |
---|
2047 | 316 FORMAT (' Sensible heatflux and momentum flux from coupled ', & |
---|
2048 | 'atmosphere model') |
---|
2049 | 317 FORMAT (//' Lateral boundaries:'/ & |
---|
2050 | ' left/right: ',A/ & |
---|
2051 | ' north/south: ',A) |
---|
2052 | 318 FORMAT (/' use_cmax: ',L1 / & |
---|
2053 | ' pt damping layer width = ',F8.2,' m, pt ', & |
---|
2054 | 'damping factor = ',F6.4) |
---|
2055 | 319 FORMAT (' turbulence recycling at inflow switched on'/ & |
---|
2056 | ' width of recycling domain: ',F7.1,' m grid index: ',I4/ & |
---|
2057 | ' inflow damping height: ',F6.1,' m width: ',F6.1,' m') |
---|
2058 | 320 FORMAT (' Predefined constant momentumflux: u: ',F9.6,' m**2/s**2'/ & |
---|
2059 | ' v: ',F9.6,' m**2/s**2') |
---|
2060 | 321 FORMAT (//' Initial profiles:'/ & |
---|
2061 | ' ----------------') |
---|
2062 | 322 FORMAT (' turbulence recycling at inflow switched on'/ & |
---|
2063 | ' y shift of the recycled inflow turbulence switched on'/ & |
---|
2064 | ' width of recycling domain: ',F7.1,' m grid index: ',I4/ & |
---|
2065 | ' inflow damping height: ',F6.1,' m width: ',F6.1,' m'/q) |
---|
2066 | 325 FORMAT (//' List output:'/ & |
---|
2067 | ' -----------'// & |
---|
2068 | ' 1D-Profiles:'/ & |
---|
2069 | ' Output every ',F8.2,' s') |
---|
2070 | 326 FORMAT (' Time averaged over ',F8.2,' s'/ & |
---|
2071 | ' Averaging input every ',F8.2,' s') |
---|
2072 | 330 FORMAT (//' Data output:'/ & |
---|
2073 | ' -----------'/) |
---|
2074 | 331 FORMAT (/' 1D-Profiles:') |
---|
2075 | 332 FORMAT (/' ',A) |
---|
2076 | 333 FORMAT (' Output every ',F8.2,' s',/ & |
---|
2077 | ' Time averaged over ',F8.2,' s'/ & |
---|
2078 | ' Averaging input every ',F8.2,' s') |
---|
2079 | 334 FORMAT (/' 2D-Arrays',A,':') |
---|
2080 | 335 FORMAT (/' ',A2,'-cross-section Arrays: ',A/ & |
---|
2081 | ' Output every ',F8.2,' s ',A/ & |
---|
2082 | ' Cross sections at ',A1,' = ',A/ & |
---|
2083 | ' scalar-coordinates: ',A,' m'/) |
---|
2084 | 336 FORMAT (/' 3D-Arrays',A,':') |
---|
2085 | 337 FORMAT (/' Arrays: ',A/ & |
---|
2086 | ' Output every ',F8.2,' s ',A/ & |
---|
2087 | ' Upper output limit at ',F8.2,' m (GP ',I4,')'/) |
---|
2088 | 339 FORMAT (' No output during initial ',F8.2,' s') |
---|
2089 | 340 FORMAT (/' Time series:') |
---|
2090 | 341 FORMAT (' Output every ',F8.2,' s'/) |
---|
2091 | 342 FORMAT (/' ',A2,'-cross-section Arrays: ',A/ & |
---|
2092 | ' Output every ',F8.2,' s ',A/ & |
---|
2093 | ' Time averaged over ',F8.2,' s'/ & |
---|
2094 | ' Averaging input every ',F8.2,' s'/ & |
---|
2095 | ' Cross sections at ',A1,' = ',A/ & |
---|
2096 | ' scalar-coordinates: ',A,' m'/) |
---|
2097 | 343 FORMAT (/' Arrays: ',A/ & |
---|
2098 | ' Output every ',F8.2,' s ',A/ & |
---|
2099 | ' Time averaged over ',F8.2,' s'/ & |
---|
2100 | ' Averaging input every ',F8.2,' s'/ & |
---|
2101 | ' Upper output limit at ',F8.2,' m (GP ',I4,')'/) |
---|
2102 | 344 FORMAT (' Output format: ',A/) |
---|
2103 | 345 FORMAT (/' Scaling lengths for output locations of all subsequent mask IDs:',/ & |
---|
2104 | ' mask_scale_x (in x-direction): ',F9.3, ' m',/ & |
---|
2105 | ' mask_scale_y (in y-direction): ',F9.3, ' m',/ & |
---|
2106 | ' mask_scale_z (in z-direction): ',F9.3, ' m' ) |
---|
2107 | 346 FORMAT (/' Masked data output',A,' for mask ID ',I2, ':') |
---|
2108 | 347 FORMAT (' Variables: ',A/ & |
---|
2109 | ' Output every ',F8.2,' s') |
---|
2110 | 348 FORMAT (' Variables: ',A/ & |
---|
2111 | ' Output every ',F8.2,' s'/ & |
---|
2112 | ' Time averaged over ',F8.2,' s'/ & |
---|
2113 | ' Averaging input every ',F8.2,' s') |
---|
2114 | 349 FORMAT (/' Output locations in ',A,'-direction in multiples of ', & |
---|
2115 | 'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ & |
---|
2116 | 13(' ',8(F8.2,',')/) ) |
---|
2117 | 350 FORMAT (/' Output locations in ',A,'-direction: ', & |
---|
2118 | 'all gridpoints along ',A,'-direction (default).' ) |
---|
2119 | 351 FORMAT (/' Output locations in ',A,'-direction in multiples of ', & |
---|
2120 | 'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ & |
---|
2121 | ' loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 ) |
---|
2122 | 352 FORMAT (/' Number of output time levels allowed: ',I3 /) |
---|
2123 | 353 FORMAT (/' Number of output time levels allowed: unlimited' /) |
---|
2124 | #if defined( __dvrp_graphics ) |
---|
2125 | 360 FORMAT (' Plot-Sequence with dvrp-software:'/ & |
---|
2126 | ' Output every ',F7.1,' s'/ & |
---|
2127 | ' Output mode: ',A/ & |
---|
2128 | ' Host / User: ',A,' / ',A/ & |
---|
2129 | ' Directory: ',A// & |
---|
2130 | ' The sequence contains:') |
---|
2131 | 361 FORMAT (/' Isosurface of "',A,'" Threshold value: ', E12.3/ & |
---|
2132 | ' Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)') |
---|
2133 | 362 FORMAT (/' Slicer plane ',A/ & |
---|
2134 | ' Slicer limits: [',F6.2,',',F6.2,']') |
---|
2135 | 363 FORMAT (/' Particles'/ & |
---|
2136 | ' particle size: ',F7.2,' m') |
---|
2137 | 364 FORMAT (' particle ',A,' controlled by "',A,'" with interval [', & |
---|
2138 | F6.2,',',F6.2,']') |
---|
2139 | 365 FORMAT (/' Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ & |
---|
2140 | ' Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, & |
---|
2141 | ')'/ & |
---|
2142 | ' Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ & |
---|
2143 | ' from y = ',F9.1,' m to y = ',F9.1,' m') |
---|
2144 | 366 FORMAT (/' Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)') |
---|
2145 | 367 FORMAT (' Polygon reduction for topography: cluster_size = ', I1) |
---|
2146 | #endif |
---|
2147 | #if defined( __spectra ) |
---|
2148 | 370 FORMAT (' Spectra:') |
---|
2149 | 371 FORMAT (' Output every ',F7.1,' s'/) |
---|
2150 | 372 FORMAT (' Arrays: ', 10(A5,',')/ & |
---|
2151 | ' Directions: ', 10(A5,',')/ & |
---|
2152 | ' height levels k = ', 20(I3,',')/ & |
---|
2153 | ' ', 20(I3,',')/ & |
---|
2154 | ' ', 20(I3,',')/ & |
---|
2155 | ' ', 20(I3,',')/ & |
---|
2156 | ' ', 19(I3,','),I3,'.'/ & |
---|
2157 | ' height levels selected for standard plot:'/ & |
---|
2158 | ' k = ', 20(I3,',')/ & |
---|
2159 | ' ', 20(I3,',')/ & |
---|
2160 | ' ', 20(I3,',')/ & |
---|
2161 | ' ', 20(I3,',')/ & |
---|
2162 | ' ', 19(I3,','),I3,'.'/ & |
---|
2163 | ' Time averaged over ', F7.1, ' s,' / & |
---|
2164 | ' Profiles for the time averaging are taken every ', & |
---|
2165 | F6.1,' s') |
---|
2166 | #endif |
---|
2167 | 400 FORMAT (//' Physical quantities:'/ & |
---|
2168 | ' -------------------'/) |
---|
2169 | 410 FORMAT (' Geograph. latitude : phi = ',F4.1,' degr'/ & |
---|
2170 | ' Angular velocity : omega = ',E9.3,' rad/s'/ & |
---|
2171 | ' Coriolis parameter : f = ',F9.6,' 1/s'/ & |
---|
2172 | ' f* = ',F9.6,' 1/s') |
---|
2173 | 411 FORMAT (/' Gravity : g = ',F4.1,' m/s**2') |
---|
2174 | 412 FORMAT (/' Reference state used in buoyancy terms: ',A) |
---|
2175 | 413 FORMAT (' Reference density in buoyancy terms: ',F8.3,' kg/m**3') |
---|
2176 | 414 FORMAT (' Reference temperature in buoyancy terms: ',F8.4,' K') |
---|
2177 | 415 FORMAT (/' Cloud physics parameters:'/ & |
---|
2178 | ' ------------------------'/) |
---|
2179 | 416 FORMAT (' Surface pressure : p_0 = ',F7.2,' hPa'/ & |
---|
2180 | ' Gas constant : R = ',F5.1,' J/(kg K)'/ & |
---|
2181 | ' Density of air : rho_0 = ',F5.3,' kg/m**3'/ & |
---|
2182 | ' Specific heat cap. : c_p = ',F6.1,' J/(kg K)'/ & |
---|
2183 | ' Vapourization heat : L_v = ',E8.2,' J/kg') |
---|
2184 | 417 FORMAT (' Geograph. longitude : lambda = ',F4.1,' degr') |
---|
2185 | 418 FORMAT (/' Day of the year at model start : day_init = ',I3 & |
---|
2186 | /' UTC time at model start : time_utc_init = ',F7.1' s') |
---|
2187 | 419 FORMAT (//' Land surface model information:'/ & |
---|
2188 | ' ------------------------------'/) |
---|
2189 | 420 FORMAT (/' Characteristic levels of the initial temperature profile:'// & |
---|
2190 | ' Height: ',A,' m'/ & |
---|
2191 | ' Temperature: ',A,' K'/ & |
---|
2192 | ' Gradient: ',A,' K/100m'/ & |
---|
2193 | ' Gridpoint: ',A) |
---|
2194 | 421 FORMAT (/' Characteristic levels of the initial humidity profile:'// & |
---|
2195 | ' Height: ',A,' m'/ & |
---|
2196 | ' Humidity: ',A,' kg/kg'/ & |
---|
2197 | ' Gradient: ',A,' (kg/kg)/100m'/ & |
---|
2198 | ' Gridpoint: ',A) |
---|
2199 | 422 FORMAT (/' Characteristic levels of the initial scalar profile:'// & |
---|
2200 | ' Height: ',A,' m'/ & |
---|
2201 | ' Scalar concentration: ',A,' kg/m**3'/ & |
---|
2202 | ' Gradient: ',A,' (kg/m**3)/100m'/ & |
---|
2203 | ' Gridpoint: ',A) |
---|
2204 | 423 FORMAT (/' Characteristic levels of the geo. wind component ug:'// & |
---|
2205 | ' Height: ',A,' m'/ & |
---|
2206 | ' ug: ',A,' m/s'/ & |
---|
2207 | ' Gradient: ',A,' 1/100s'/ & |
---|
2208 | ' Gridpoint: ',A) |
---|
2209 | 424 FORMAT (/' Characteristic levels of the geo. wind component vg:'// & |
---|
2210 | ' Height: ',A,' m'/ & |
---|
2211 | ' vg: ',A,' m/s'/ & |
---|
2212 | ' Gradient: ',A,' 1/100s'/ & |
---|
2213 | ' Gridpoint: ',A) |
---|
2214 | 425 FORMAT (/' Characteristic levels of the initial salinity profile:'// & |
---|
2215 | ' Height: ',A,' m'/ & |
---|
2216 | ' Salinity: ',A,' psu'/ & |
---|
2217 | ' Gradient: ',A,' psu/100m'/ & |
---|
2218 | ' Gridpoint: ',A) |
---|
2219 | 426 FORMAT (/' Characteristic levels of the subsidence/ascent profile:'// & |
---|
2220 | ' Height: ',A,' m'/ & |
---|
2221 | ' w_subs: ',A,' m/s'/ & |
---|
2222 | ' Gradient: ',A,' (m/s)/100m'/ & |
---|
2223 | ' Gridpoint: ',A) |
---|
2224 | 427 FORMAT (/' Initial wind profiles (u,v) are interpolated from given'// & |
---|
2225 | ' profiles') |
---|
2226 | 428 FORMAT (/' Initial profiles (u, v, pt, q) are taken from file '/ & |
---|
2227 | ' NUDGING_DATA') |
---|
2228 | 430 FORMAT (//' Cloud physics quantities / methods:'/ & |
---|
2229 | ' ----------------------------------'/) |
---|
2230 | 431 FORMAT (' Humidity is treated as purely passive scalar (no condensati', & |
---|
2231 | 'on)') |
---|
2232 | 432 FORMAT (' Bulk scheme with liquid water potential temperature and'/ & |
---|
2233 | ' total water content is used.'/ & |
---|
2234 | ' Condensation is parameterized via 0% - or 100% scheme.') |
---|
2235 | 433 FORMAT (' Cloud droplets treated explicitly using the Lagrangian part', & |
---|
2236 | 'icle model') |
---|
2237 | 434 FORMAT (' Curvature and solution effecs are considered for growth of', & |
---|
2238 | ' droplets < 1.0E-6 m') |
---|
2239 | 435 FORMAT (' Droplet collision is handled by ',A,'-kernel') |
---|
2240 | 436 FORMAT (' Fast kernel with fixed radius- and dissipation classes ', & |
---|
2241 | 'are used'/ & |
---|
2242 | ' number of radius classes: ',I3,' interval ', & |
---|
2243 | '[1.0E-6,2.0E-4] m'/ & |
---|
2244 | ' number of dissipation classes: ',I2,' interval ', & |
---|
2245 | '[0,1000] cm**2/s**3') |
---|
2246 | 437 FORMAT (' Droplet collision is switched off') |
---|
2247 | 438 FORMAT (' --> Land surface type : ',A,/ & |
---|
2248 | ' --> Soil porosity type : ',A) |
---|
2249 | 439 FORMAT (/' Initial soil temperature and moisture profile:'// & |
---|
2250 | ' Height: ',A,' m'/ & |
---|
2251 | ' Temperature: ',A,' K'/ & |
---|
2252 | ' Moisture: ',A,' m**3/m**3'/ & |
---|
2253 | ' Root fraction: ',A,' '/ & |
---|
2254 | ' Gridpoint: ',A) |
---|
2255 | 440 FORMAT (/' --> Dewfall is allowed (default)') |
---|
2256 | 441 FORMAT (' --> Dewfall is inhibited') |
---|
2257 | 442 FORMAT (' --> Soil bottom is closed (water content is conserved, default)') |
---|
2258 | 443 FORMAT (' --> Soil bottom is open (water content is not conserved)') |
---|
2259 | 444 FORMAT (//' Radiation model information:'/ & |
---|
2260 | ' ----------------------------'/) |
---|
2261 | 445 FORMAT (' --> Using constant net radiation: net_radiation = ', F6.2, ' W/m**2') |
---|
2262 | 446 FORMAT (' --> Simple radiation scheme for clear sky is used (no clouds,', & |
---|
2263 | ' default)') |
---|
2264 | 447 FORMAT (' --> Radiation scheme:', A) |
---|
2265 | 448 FORMAT (/' Surface albedo: albedo = ', F5.3) |
---|
2266 | 449 FORMAT (' Timestep: dt_radiation = ', F5.2, ' s') |
---|
2267 | |
---|
2268 | 450 FORMAT (//' LES / Turbulence quantities:'/ & |
---|
2269 | ' ---------------------------'/) |
---|
2270 | 451 FORMAT (' Diffusion coefficients are constant:'/ & |
---|
2271 | ' Km = ',F6.2,' m**2/s Kh = ',F6.2,' m**2/s Pr = ',F5.2) |
---|
2272 | 453 FORMAT (' Mixing length is limited to ',F4.2,' * z') |
---|
2273 | 454 FORMAT (' TKE is not allowed to fall below ',E9.2,' (m/s)**2') |
---|
2274 | 455 FORMAT (' initial TKE is prescribed as ',E9.2,' (m/s)**2') |
---|
2275 | 470 FORMAT (//' Actions during the simulation:'/ & |
---|
2276 | ' -----------------------------'/) |
---|
2277 | 471 FORMAT (' Disturbance impulse (u,v) every : ',F6.2,' s'/ & |
---|
2278 | ' Disturbance amplitude : ',F4.2, ' m/s'/ & |
---|
2279 | ' Lower disturbance level : ',F8.2,' m (GP ',I4,')'/ & |
---|
2280 | ' Upper disturbance level : ',F8.2,' m (GP ',I4,')') |
---|
2281 | 472 FORMAT (' Disturbances continued during the run from i/j =',I4, & |
---|
2282 | ' to i/j =',I4) |
---|
2283 | 473 FORMAT (' Disturbances cease as soon as the disturbance energy exceeds',& |
---|
2284 | 1X,F5.3, ' m**2/s**2') |
---|
2285 | 474 FORMAT (' Random number generator used : ',A/) |
---|
2286 | 475 FORMAT (' The surface temperature is increased (or decreased, ', & |
---|
2287 | 'respectively, if'/ & |
---|
2288 | ' the value is negative) by ',F5.2,' K at the beginning of the',& |
---|
2289 | ' 3D-simulation'/) |
---|
2290 | 476 FORMAT (' The surface humidity is increased (or decreased, ',& |
---|
2291 | 'respectively, if the'/ & |
---|
2292 | ' value is negative) by ',E8.1,' kg/kg at the beginning of', & |
---|
2293 | ' the 3D-simulation'/) |
---|
2294 | 477 FORMAT (' The scalar value is increased at the surface (or decreased, ',& |
---|
2295 | 'respectively, if the'/ & |
---|
2296 | ' value is negative) by ',E8.1,' kg/m**3 at the beginning of', & |
---|
2297 | ' the 3D-simulation'/) |
---|
2298 | 480 FORMAT (' Particles:'/ & |
---|
2299 | ' ---------'// & |
---|
2300 | ' Particle advection is active (switched on at t = ', F7.1, & |
---|
2301 | ' s)'/ & |
---|
2302 | ' Start of new particle generations every ',F6.1,' s'/ & |
---|
2303 | ' Boundary conditions: left/right: ', A, ' north/south: ', A/& |
---|
2304 | ' bottom: ', A, ' top: ', A/& |
---|
2305 | ' Maximum particle age: ',F9.1,' s'/ & |
---|
2306 | ' Advection stopped at t = ',F9.1,' s'/) |
---|
2307 | 481 FORMAT (' Particles have random start positions'/) |
---|
2308 | 482 FORMAT (' Particles are advected only horizontally'/) |
---|
2309 | 483 FORMAT (' Particles have tails with a maximum of ',I3,' points') |
---|
2310 | 484 FORMAT (' Number of tails of the total domain: ',I10/ & |
---|
2311 | ' Minimum distance between tailpoints: ',F8.2,' m'/ & |
---|
2312 | ' Maximum age of the end of the tail: ',F8.2,' s') |
---|
2313 | 485 FORMAT (' Particle data are written on file every ', F9.1, ' s') |
---|
2314 | 486 FORMAT (' Particle statistics are written on file'/) |
---|
2315 | 487 FORMAT (' Number of particle groups: ',I2/) |
---|
2316 | 488 FORMAT (' SGS velocity components are used for particle advection'/ & |
---|
2317 | ' minimum timestep for advection: ', F7.5/) |
---|
2318 | 489 FORMAT (' Number of particles simultaneously released at each ', & |
---|
2319 | 'point: ', I5/) |
---|
2320 | 490 FORMAT (' Particle group ',I2,':'/ & |
---|
2321 | ' Particle radius: ',E10.3, 'm') |
---|
2322 | 491 FORMAT (' Particle inertia is activated'/ & |
---|
2323 | ' density_ratio (rho_fluid/rho_particle) = ',F5.3/) |
---|
2324 | 492 FORMAT (' Particles are advected only passively (no inertia)'/) |
---|
2325 | 493 FORMAT (' Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/& |
---|
2326 | ' y:',F8.1,' - ',F8.1,' m'/& |
---|
2327 | ' z:',F8.1,' - ',F8.1,' m'/& |
---|
2328 | ' Particle distances: dx = ',F8.1,' m dy = ',F8.1, & |
---|
2329 | ' m dz = ',F8.1,' m'/) |
---|
2330 | 494 FORMAT (' Output of particle time series in NetCDF format every ', & |
---|
2331 | F8.2,' s'/) |
---|
2332 | 495 FORMAT (' Number of particles in total domain: ',I10/) |
---|
2333 | 496 FORMAT (' Initial vertical particle positions are interpreted ', & |
---|
2334 | 'as relative to the given topography') |
---|
2335 | 500 FORMAT (//' 1D-Model parameters:'/ & |
---|
2336 | ' -------------------'// & |
---|
2337 | ' Simulation time: ',F8.1,' s'/ & |
---|
2338 | ' Run-controll output every: ',F8.1,' s'/ & |
---|
2339 | ' Vertical profile output every: ',F8.1,' s'/ & |
---|
2340 | ' Mixing length calculation: ',A/ & |
---|
2341 | ' Dissipation calculation: ',A/) |
---|
2342 | 502 FORMAT (' Damping layer starts from ',F7.1,' m (GP ',I4,')'/) |
---|
2343 | 503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order') |
---|
2344 | 504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order') |
---|
2345 | 505 FORMAT (' Precipitation parameterization via Seifert-Beheng-Scheme') |
---|
2346 | 506 FORMAT (' Drizzle parameterization via Stokes law') |
---|
2347 | 507 FORMAT (' Turbulence effects on precipitation process') |
---|
2348 | 508 FORMAT (' Ventilation effects on evaporation of rain drops') |
---|
2349 | 509 FORMAT (' Slope limiter used for sedimentation process') |
---|
2350 | 510 FORMAT (' Droplet density : N_c = ',F6.1,' 1/cm**3') |
---|
2351 | 511 FORMAT (' Sedimentation Courant number: '/& |
---|
2352 | ' C_s = ',F3.1,' ') |
---|
2353 | 512 FORMAT (/' Date: ',A8,6X,'Run: ',A20/ & |
---|
2354 | ' Time: ',A8,6X,'Run-No.: ',I2.2/ & |
---|
2355 | ' Run on host: ',A10,6X,'En-No.: ',I2.2) |
---|
2356 | 513 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order ' // & |
---|
2357 | '+ monotonic adjustment') |
---|
2358 | |
---|
2359 | |
---|
2360 | END SUBROUTINE header |
---|