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

Last change on this file since 2014 was 2001, checked in by knoop, 8 years ago

last commit documented

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