source: palm/trunk/SOURCE/header.f90 @ 1817

Last change on this file since 1817 was 1817, checked in by maronga, 8 years ago

further modularization of land surface model

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