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

Last change on this file since 2050 was 2050, checked in by gronemeier, 7 years ago

Implement turbulent outflow condition

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