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

Last change on this file since 1359 was 1359, checked in by hoffmann, 10 years ago

new Lagrangian particle structure integrated

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