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

Last change on this file since 2073 was 2073, checked in by raasch, 7 years ago

mostly openmp related bugfixes

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