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

Last change on this file since 978 was 978, checked in by fricke, 12 years ago

merge fricke branch back into trunk

  • Property svn:keywords set to Id
File size: 79.7 KB
Line 
1 SUBROUTINE header
2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
6! -km_damp_max, outflow_damping_width
7! +pt_damping_factor, pt_damping_width
8! +z0h
9!
10! Former revisions:
11! -----------------
12! $Id: header.f90 978 2012-08-09 08:28:32Z fricke $
13!
14! 964 2012-07-26 09:14:24Z raasch
15! output of profil-related quantities removed
16!
17! 940 2012-07-09 14:31:00Z raasch
18! Output in case of simulations for pure neutral stratification (no pt-equation
19! solved)
20!
21! 927 2012-06-06 19:15:04Z raasch
22! output of masking_method for mg-solver
23!
24! 868 2012-03-28 12:21:07Z raasch
25! translation velocity in Galilean transformation changed to 0.6 * ug
26!
27! 833 2012-02-22 08:55:55Z maronga
28! Adjusted format for leaf area density
29!
30! 828 2012-02-21 12:00:36Z raasch
31! output of dissipation_classes + radius_classes
32!
33! 825 2012-02-19 03:03:44Z raasch
34! Output of cloud physics parameters/quantities complemented and restructured
35!
36! 767 2011-10-14 06:39:12Z raasch
37! Output of given initial u,v-profiles
38!
39! 759 2011-09-15 13:58:31Z raasch
40! output of maximum number of parallel io streams
41!
42! 707 2011-03-29 11:39:40Z raasch
43! bc_lr/ns replaced by bc_lr/ns_cyc
44!
45! 667 2010-12-23 12:06:00Z suehring/gryschka
46! Output of advection scheme.
47! Modified output of Prandtl-layer height.
48!
49! 580 2010-10-05 13:59:11Z heinze
50! Renaming of ws_vertical_gradient to subs_vertical_gradient,
51! ws_vertical_gradient_level to subs_vertical_gradient_level and
52! ws_vertical_gradient_level_ind to subs_vertical_gradient_level_i
53!
54! 493 2010-03-01 08:30:24Z raasch
55! NetCDF data output format extendend for NetCDF4/HDF5
56!
57! 449 2010-02-02 11:23:59Z raasch
58! +large scale vertical motion (subsidence/ascent)
59! Bugfix: index problem concerning gradient_level indices removed
60!
61! 410 2009-12-04 17:05:40Z letzel
62! Masked data output: + dt_domask, mask_01~20_x|y|z, mask_01~20_x|y|z_loop,
63! mask_scale|_x|y|z, masks, skip_time_domask
64!
65! 346 2009-07-06 10:13:41Z raasch
66! initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill'
67! Coupling with independent precursor runs.
68! Output of messages replaced by message handling routine.
69! Output of several additional dvr parameters
70! +canyon_height, canyon_width_x, canyon_width_y, canyon_wall_left,
71! canyon_wall_south, conserve_volume_flow_mode, dp_external, dp_level_b,
72! dp_smooth, dpdxy, u_bulk, v_bulk
73! topography_grid_convention moved from user_header
74! small bugfix concerning 3d 64bit netcdf output format
75!
76! 206 2008-10-13 14:59:11Z raasch
77! Bugfix: error in zu index in case of section_xy = -1
78!
79! 198 2008-09-17 08:55:28Z raasch
80! Format adjustments allowing output of larger revision numbers
81!
82! 197 2008-09-16 15:29:03Z raasch
83! allow 100 spectra levels instead of 10 for consistency with
84! define_netcdf_header,
85! bugfix in the output of the characteristic levels of potential temperature,
86! geostrophic wind, scalar concentration, humidity and leaf area density,
87! output of turbulence recycling informations
88!
89! 138 2007-11-28 10:03:58Z letzel
90! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
91! Allow two instead of one digit to specify isosurface and slicer variables.
92! Output of sorting frequency of particles
93!
94! 108 2007-08-24 15:10:38Z letzel
95! Output of informations for coupled model runs (boundary conditions etc.)
96! + output of momentumfluxes at the top boundary
97! Rayleigh damping for ocean, e_init
98!
99! 97 2007-06-21 08:23:15Z raasch
100! Adjustments for the ocean version.
101! use_pt_reference renamed use_reference
102!
103! 87 2007-05-22 15:46:47Z raasch
104! Bugfix: output of use_upstream_for_tke
105!
106! 82 2007-04-16 15:40:52Z raasch
107! Preprocessor strings for different linux clusters changed to "lc",
108! routine local_flush is used for buffer flushing
109!
110! 76 2007-03-29 00:58:32Z raasch
111! Output of netcdf_64bit_3d, particles-package is now part of the default code,
112! output of the loop optimization method, moisture renamed humidity,
113! output of subversion revision number
114!
115! 19 2007-02-23 04:53:48Z raasch
116! Output of scalar flux applied at top boundary
117!
118! RCS Log replace by Id keyword, revision history cleaned up
119!
120! Revision 1.63  2006/08/22 13:53:13  raasch
121! Output of dz_max
122!
123! Revision 1.1  1997/08/11 06:17:20  raasch
124! Initial revision
125!
126!
127! Description:
128! ------------
129! Writing a header with all important informations about the actual run.
130! This subroutine is called three times, two times at the beginning
131! (writing information on files RUN_CONTROL and HEADER) and one time at the
132! end of the run, then writing additional information about CPU-usage on file
133! header.
134!-----------------------------------------------------------------------------!
135
136    USE arrays_3d
137    USE control_parameters
138    USE cloud_parameters
139    USE cpulog
140    USE dvrp_variables
141    USE grid_variables
142    USE indices
143    USE model_1d
144    USE particle_attributes
145    USE pegrid
146    USE subsidence_mod
147    USE spectrum
148
149    IMPLICIT NONE
150
151    CHARACTER (LEN=1)  ::  prec
152    CHARACTER (LEN=2)  ::  do2d_mode
153    CHARACTER (LEN=5)  ::  section_chr
154    CHARACTER (LEN=9)  ::  time_to_string
155    CHARACTER (LEN=10) ::  coor_chr, host_chr
156    CHARACTER (LEN=16) ::  begin_chr
157    CHARACTER (LEN=23) ::  ver_rev
158    CHARACTER (LEN=40) ::  output_format
159    CHARACTER (LEN=70) ::  char1, char2, dopr_chr, &
160                           do2d_xy, do2d_xz, do2d_yz, do3d_chr, &
161                           domask_chr, run_classification
162    CHARACTER (LEN=86) ::  coordinates, gradients, learde, slices,  &
163                           temperatures, ugcomponent, vgcomponent
164    CHARACTER (LEN=85) ::  roben, runten
165
166    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)
167
168    INTEGER ::  av, bh, blx, bly, bxl, bxr, byn, bys, ch, count, cwx, cwy,  &
169         cxl, cxr, cyn, cys, dim, i, ihost, io, j, l, ll, m, mpi_type
170    REAL    ::  cpuseconds_per_simulated_second
171
172!
173!-- Open the output file. At the end of the simulation, output is directed
174!-- to unit 19.
175    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
176         .NOT. simulated_time_at_begin /= simulated_time )  THEN
177       io = 15   !  header output on file RUN_CONTROL
178    ELSE
179       io = 19   !  header output on file HEADER
180    ENDIF
181    CALL check_open( io )
182
183!
184!-- At the end of the run, output file (HEADER) will be rewritten with
185!-- new informations
186    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
187
188!
189!-- Determine kind of model run
190    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
191       run_classification = '3D - restart run'
192    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
193       run_classification = '3D - run with cyclic fill of 3D - prerun data'
194    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
195       run_classification = '3D - run without 1D - prerun'
196    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
197       run_classification = '3D - run with 1D - prerun'
198    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
199       run_classification = '3D - run initialized by user'
200    ELSE
201       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
202       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
203    ENDIF
204    IF ( ocean )  THEN
205       run_classification = 'ocean - ' // run_classification
206    ELSE
207       run_classification = 'atmosphere - ' // run_classification
208    ENDIF
209
210!
211!-- Run-identification, date, time, host
212    host_chr = host(1:10)
213    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
214    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
215    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
216#if defined( __mpi2 )
217       mpi_type = 2
218#else
219       mpi_type = 1
220#endif
221       WRITE ( io, 101 )  mpi_type, coupling_mode
222    ENDIF
223    WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr, &
224                       ADJUSTR( host_chr )
225#if defined( __parallel )
226    IF ( npex == -1  .AND.  pdims(2) /= 1 )  THEN
227       char1 = 'calculated'
228    ELSEIF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.  &
229               host(1:2) == 'lc' )  .AND.                          &
230             npex == -1  .AND.  pdims(2) == 1 )  THEN
231       char1 = 'forced'
232    ELSE
233       char1 = 'predefined'
234    ENDIF
235    IF ( threads_per_task == 1 )  THEN
236       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
237    ELSE
238       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
239                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
240    ENDIF
241    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
242           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
243         npex == -1  .AND.  pdims(2) == 1 )                      &
244    THEN
245       WRITE ( io, 106 )
246    ELSEIF ( pdims(2) == 1 )  THEN
247       WRITE ( io, 107 )  'x'
248    ELSEIF ( pdims(1) == 1 )  THEN
249       WRITE ( io, 107 )  'y'
250    ENDIF
251    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
252    IF ( numprocs /= maximum_parallel_io_streams )  THEN
253       WRITE ( io, 108 )  maximum_parallel_io_streams
254    ENDIF
255#endif
256    WRITE ( io, 99 )
257
258!
259!-- Numerical schemes
260    WRITE ( io, 110 )
261    IF ( psolver(1:7) == 'poisfft' )  THEN
262       WRITE ( io, 111 )  TRIM( fft_method )
263       IF ( psolver == 'poisfft_hybrid' )  WRITE ( io, 138 )
264    ELSEIF ( psolver == 'sor' )  THEN
265       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
266    ELSEIF ( psolver == 'multigrid' )  THEN
267       WRITE ( io, 135 )  cycle_mg, maximum_grid_level, ngsrb
268       IF ( mg_cycles == -1 )  THEN
269          WRITE ( io, 140 )  residual_limit
270       ELSE
271          WRITE ( io, 141 )  mg_cycles
272       ENDIF
273       IF ( mg_switch_to_pe0_level == 0 )  THEN
274          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
275                             nzt_mg(1)
276       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
277          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
278                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
279                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
280                             nzt_mg(mg_switch_to_pe0_level),    &
281                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
282                             nzt_mg(1)
283       ENDIF
284       IF ( masking_method )  WRITE ( io, 144 )
285    ENDIF
286    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
287    THEN
288       WRITE ( io, 142 )
289    ENDIF
290
291    IF ( momentum_advec == 'pw-scheme' )  THEN
292       WRITE ( io, 113 )
293    ELSEIF (momentum_advec == 'ws-scheme' ) THEN
294       WRITE ( io, 503 )
295    ELSEIF (momentum_advec == 'ups-scheme' ) THEN
296       WRITE ( io, 114 )
297       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
298       IF ( overshoot_limit_u /= 0.0  .OR.  overshoot_limit_v /= 0.0  .OR. &
299            overshoot_limit_w /= 0.0 )  THEN
300          WRITE ( io, 127 )  overshoot_limit_u, overshoot_limit_v, &
301                             overshoot_limit_w
302       ENDIF
303       IF ( ups_limit_u /= 0.0  .OR.  ups_limit_v /= 0.0  .OR. &
304            ups_limit_w /= 0.0 )                               &
305       THEN
306          WRITE ( io, 125 )  ups_limit_u, ups_limit_v, ups_limit_w
307       ENDIF
308       IF ( long_filter_factor /= 0.0 )  WRITE ( io, 115 )  long_filter_factor
309    ENDIF
310    IF ( scalar_advec == 'pw-scheme' )  THEN
311       WRITE ( io, 116 )
312    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
313       WRITE ( io, 504 )
314    ELSEIF ( scalar_advec == 'ups-scheme' )  THEN
315       WRITE ( io, 117 )
316       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
317       IF ( overshoot_limit_e /= 0.0  .OR.  overshoot_limit_pt /= 0.0 )  THEN
318          WRITE ( io, 128 )  overshoot_limit_e, overshoot_limit_pt
319       ENDIF
320       IF ( ups_limit_e /= 0.0  .OR.  ups_limit_pt /= 0.0 )  THEN
321          WRITE ( io, 126 )  ups_limit_e, ups_limit_pt
322       ENDIF
323    ELSE
324       WRITE ( io, 118 )
325    ENDIF
326
327    WRITE ( io, 139 )  TRIM( loop_optimization )
328
329    IF ( galilei_transformation )  THEN
330       IF ( use_ug_for_galilei_tr )  THEN
331          char1 = '0.6 * geostrophic wind'
332       ELSE
333          char1 = 'mean wind in model domain'
334       ENDIF
335       IF ( simulated_time_at_begin == simulated_time )  THEN
336          char2 = 'at the start of the run'
337       ELSE
338          char2 = 'at the end of the run'
339       ENDIF
340       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ), &
341                          advected_distance_x/1000.0, advected_distance_y/1000.0
342    ENDIF
343    IF ( timestep_scheme == 'leapfrog' )  THEN
344       WRITE ( io, 120 )
345    ELSEIF ( timestep_scheme == 'leapfrog+euler' )  THEN
346       WRITE ( io, 121 )
347    ELSE
348       WRITE ( io, 122 )  timestep_scheme
349    ENDIF
350    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
351    IF ( rayleigh_damping_factor /= 0.0 )  THEN
352       IF ( .NOT. ocean )  THEN
353          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
354               rayleigh_damping_factor
355       ELSE
356          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
357               rayleigh_damping_factor
358       ENDIF
359    ENDIF
360    IF ( neutral )  WRITE ( io, 131 )  pt_surface
361    IF ( humidity )  THEN
362       IF ( .NOT. cloud_physics )  THEN
363          WRITE ( io, 129 )
364       ELSE
365          WRITE ( io, 130 )
366       ENDIF
367    ENDIF
368    IF ( passive_scalar )  WRITE ( io, 134 )
369    IF ( conserve_volume_flow )  THEN
370       WRITE ( io, 150 )  conserve_volume_flow_mode
371       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
372          WRITE ( io, 151 )  u_bulk, v_bulk
373       ENDIF
374    ELSEIF ( dp_external )  THEN
375       IF ( dp_smooth )  THEN
376          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
377       ELSE
378          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
379       ENDIF
380    ENDIF
381    IF ( large_scale_subsidence )  THEN
382        WRITE ( io, 153 )
383        WRITE ( io, 154 )
384    ENDIF
385    WRITE ( io, 99 )
386
387!
388!-- Runtime and timestep informations
389    WRITE ( io, 200 )
390    IF ( .NOT. dt_fixed )  THEN
391       WRITE ( io, 201 )  dt_max, cfl_factor
392    ELSE
393       WRITE ( io, 202 )  dt
394    ENDIF
395    WRITE ( io, 203 )  simulated_time_at_begin, end_time
396
397    IF ( time_restart /= 9999999.9  .AND. &
398         simulated_time_at_begin == simulated_time )  THEN
399       IF ( dt_restart == 9999999.9 )  THEN
400          WRITE ( io, 204 )  ' Restart at:       ',time_restart
401       ELSE
402          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
403       ENDIF
404    ENDIF
405
406    IF ( simulated_time_at_begin /= simulated_time )  THEN
407       i = MAX ( log_point_s(10)%counts, 1 )
408       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0 )  THEN
409          cpuseconds_per_simulated_second = 0.0
410       ELSE
411          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
412                                            ( simulated_time -    &
413                                              simulated_time_at_begin )
414       ENDIF
415       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum, &
416                          log_point_s(10)%sum / REAL( i ),     &
417                          cpuseconds_per_simulated_second
418       IF ( time_restart /= 9999999.9  .AND.  time_restart < end_time )  THEN
419          IF ( dt_restart == 9999999.9 )  THEN
420             WRITE ( io, 204 )  ' Next restart at:  ',time_restart
421          ELSE
422             WRITE ( io, 205 )  ' Next restart at:  ',time_restart, dt_restart
423          ENDIF
424       ENDIF
425    ENDIF
426
427!
428!-- Start time for coupled runs, if independent precursor runs for atmosphere
429!-- and ocean are used. In this case, coupling_start_time defines the time
430!-- when the coupling is switched on.
431    IF ( coupling_start_time /= 0.0 )  THEN
432       IF ( coupling_start_time >= simulated_time_at_begin )  THEN
433          char1 = 'Precursor run for a coupled atmosphere-ocean run'
434       ELSE
435          char1 = 'Coupled atmosphere-ocean run following independent ' // &
436                  'precursor runs'
437       ENDIF
438       WRITE ( io, 207 )  char1, coupling_start_time
439    ENDIF
440
441!
442!-- Computational grid
443    IF ( .NOT. ocean )  THEN
444       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
445       IF ( dz_stretch_level_index < nzt+1 )  THEN
446          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
447                             dz_stretch_factor, dz_max
448       ENDIF
449    ELSE
450       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
451       IF ( dz_stretch_level_index > 0 )  THEN
452          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
453                             dz_stretch_factor, dz_max
454       ENDIF
455    ENDIF
456    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
457                       MIN( nnz+2, nzt+2 )
458    IF ( numprocs > 1 )  THEN
459       IF ( nxa == nx  .AND.  nya == ny  .AND.  nza == nz )  THEN
460          WRITE ( io, 255 )
461       ELSE
462          WRITE ( io, 256 )  nnx-(nxa-nx), nny-(nya-ny), nzt+2
463       ENDIF
464    ENDIF
465    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
466
467!
468!-- Topography
469    WRITE ( io, 270 )  topography
470    SELECT CASE ( TRIM( topography ) )
471
472       CASE ( 'flat' )
473          ! no actions necessary
474
475       CASE ( 'single_building' )
476          blx = INT( building_length_x / dx )
477          bly = INT( building_length_y / dy )
478          bh  = INT( building_height / dz )
479
480          IF ( building_wall_left == 9999999.9 )  THEN
481             building_wall_left = ( nx + 1 - blx ) / 2 * dx
482          ENDIF
483          bxl = INT ( building_wall_left / dx + 0.5 )
484          bxr = bxl + blx
485
486          IF ( building_wall_south == 9999999.9 )  THEN
487             building_wall_south = ( ny + 1 - bly ) / 2 * dy
488          ENDIF
489          bys = INT ( building_wall_south / dy + 0.5 )
490          byn = bys + bly
491
492          WRITE ( io, 271 )  building_length_x, building_length_y, &
493                             building_height, bxl, bxr, bys, byn
494
495       CASE ( 'single_street_canyon' )
496          ch  = NINT( canyon_height / dz )
497          IF ( canyon_width_x /= 9999999.9 )  THEN
498!
499!--          Street canyon in y direction
500             cwx = NINT( canyon_width_x / dx )
501             IF ( canyon_wall_left == 9999999.9 )  THEN
502                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
503             ENDIF
504             cxl = NINT( canyon_wall_left / dx )
505             cxr = cxl + cwx
506             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
507
508          ELSEIF ( canyon_width_y /= 9999999.9 )  THEN
509!
510!--          Street canyon in x direction
511             cwy = NINT( canyon_width_y / dy )
512             IF ( canyon_wall_south == 9999999.9 )  THEN
513                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
514             ENDIF
515             cys = NINT( canyon_wall_south / dy )
516             cyn = cys + cwy
517             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
518          ENDIF
519
520    END SELECT
521
522    IF ( TRIM( topography ) /= 'flat' )  THEN
523       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
524          IF ( TRIM( topography ) == 'single_building' .OR.  &
525               TRIM( topography ) == 'single_street_canyon' )  THEN
526             WRITE ( io, 278 )
527          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
528             WRITE ( io, 279 )
529          ENDIF
530       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
531          WRITE ( io, 278 )
532       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
533          WRITE ( io, 279 )
534       ENDIF
535    ENDIF
536
537    IF ( plant_canopy ) THEN
538
539       WRITE ( io, 280 ) canopy_mode, pch_index, drag_coefficient
540       IF ( passive_scalar ) THEN
541          WRITE ( io, 281 ) scalar_exchange_coefficient,   &
542                            leaf_surface_concentration
543       ENDIF
544
545!
546!--    Heat flux at the top of vegetation
547       WRITE ( io, 282 ) cthf
548
549!
550!--    Leaf area density profile
551!--    Building output strings, starting with surface value
552       WRITE ( learde, '(F6.4)' )  lad_surface
553       gradients = '------'
554       slices = '     0'
555       coordinates = '   0.0'
556       i = 1
557       DO  WHILE ( lad_vertical_gradient_level_ind(i) /= -9999 )
558
559          WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
560          learde = TRIM( learde ) // ' ' // TRIM( coor_chr )
561
562          WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
563          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
564
565          WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
566          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
567
568          WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
569          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
570
571          i = i + 1
572       ENDDO
573
574       WRITE ( io, 283 )  TRIM( coordinates ), TRIM( learde ), &
575                          TRIM( gradients ), TRIM( slices )
576
577    ENDIF
578
579!
580!-- Boundary conditions
581    IF ( ibc_p_b == 0 )  THEN
582       runten = 'p(0)     = 0      |'
583    ELSEIF ( ibc_p_b == 1 )  THEN
584       runten = 'p(0)     = p(1)   |'
585    ELSE
586       runten = 'p(0)     = p(1) +R|'
587    ENDIF
588    IF ( ibc_p_t == 0 )  THEN
589       roben  = 'p(nzt+1) = 0      |'
590    ELSE
591       roben  = 'p(nzt+1) = p(nzt) |'
592    ENDIF
593
594    IF ( ibc_uv_b == 0 )  THEN
595       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
596    ELSE
597       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
598    ENDIF
599    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
600       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
601    ELSEIF ( ibc_uv_t == 0 )  THEN
602       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
603    ELSE
604       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
605    ENDIF
606
607    IF ( ibc_pt_b == 0 )  THEN
608       runten = TRIM( runten ) // ' pt(0)   = pt_surface'
609    ELSEIF ( ibc_pt_b == 1 )  THEN
610       runten = TRIM( runten ) // ' pt(0)   = pt(1)'
611    ELSEIF ( ibc_pt_b == 2 )  THEN
612       runten = TRIM( runten ) // ' pt(0) = from coupled model'
613    ENDIF
614    IF ( ibc_pt_t == 0 )  THEN
615       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
616    ELSEIF( ibc_pt_t == 1 )  THEN
617       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
618    ELSEIF( ibc_pt_t == 2 )  THEN
619       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
620
621    ENDIF
622
623    WRITE ( io, 300 )  runten, roben
624
625    IF ( .NOT. constant_diffusion )  THEN
626       IF ( ibc_e_b == 1 )  THEN
627          runten = 'e(0)     = e(1)'
628       ELSE
629          runten = 'e(0)     = e(1) = (u*/0.1)**2'
630       ENDIF
631       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
632
633       WRITE ( io, 301 )  'e', runten, roben       
634
635    ENDIF
636
637    IF ( ocean )  THEN
638       runten = 'sa(0)    = sa(1)'
639       IF ( ibc_sa_t == 0 )  THEN
640          roben =  'sa(nzt+1) = sa_surface'
641       ELSE
642          roben =  'sa(nzt+1) = sa(nzt)'
643       ENDIF
644       WRITE ( io, 301 ) 'sa', runten, roben
645    ENDIF
646
647    IF ( humidity )  THEN
648       IF ( ibc_q_b == 0 )  THEN
649          runten = 'q(0)     = q_surface'
650       ELSE
651          runten = 'q(0)     = q(1)'
652       ENDIF
653       IF ( ibc_q_t == 0 )  THEN
654          roben =  'q(nzt)   = q_top'
655       ELSE
656          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
657       ENDIF
658       WRITE ( io, 301 ) 'q', runten, roben
659    ENDIF
660
661    IF ( passive_scalar )  THEN
662       IF ( ibc_q_b == 0 )  THEN
663          runten = 's(0)     = s_surface'
664       ELSE
665          runten = 's(0)     = s(1)'
666       ENDIF
667       IF ( ibc_q_t == 0 )  THEN
668          roben =  's(nzt)   = s_top'
669       ELSE
670          roben =  's(nzt)   = s(nzt-1) + ds/dz'
671       ENDIF
672       WRITE ( io, 301 ) 's', runten, roben
673    ENDIF
674
675    IF ( use_surface_fluxes )  THEN
676       WRITE ( io, 303 )
677       IF ( constant_heatflux )  THEN
678          WRITE ( io, 306 )  surface_heatflux
679          IF ( random_heatflux )  WRITE ( io, 307 )
680       ENDIF
681       IF ( humidity  .AND.  constant_waterflux )  THEN
682          WRITE ( io, 311 ) surface_waterflux
683       ENDIF
684       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
685          WRITE ( io, 313 ) surface_waterflux
686       ENDIF
687    ENDIF
688
689    IF ( use_top_fluxes )  THEN
690       WRITE ( io, 304 )
691       IF ( coupling_mode == 'uncoupled' )  THEN
692          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
693          IF ( constant_top_heatflux )  THEN
694             WRITE ( io, 306 )  top_heatflux
695          ENDIF
696       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
697          WRITE ( io, 316 )
698       ENDIF
699       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
700          WRITE ( io, 309 )  top_salinityflux
701       ENDIF
702       IF ( humidity  .OR.  passive_scalar )  THEN
703          WRITE ( io, 315 )
704       ENDIF
705    ENDIF
706
707    IF ( prandtl_layer )  THEN
708       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length, &
709                          z0h_factor*roughness_length, kappa, &
710                          rif_min, rif_max
711       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
712       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
713          WRITE ( io, 312 )
714       ENDIF
715       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
716          WRITE ( io, 314 )
717       ENDIF
718    ELSE
719       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
720          WRITE ( io, 310 )  rif_min, rif_max
721       ENDIF
722    ENDIF
723
724    WRITE ( io, 317 )  bc_lr, bc_ns
725    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
726       WRITE ( io, 318 )  pt_damping_width, pt_damping_factor       
727       IF ( turbulent_inflow )  THEN
728          WRITE ( io, 319 )  recycling_width, recycling_plane, &
729                             inflow_damping_height, inflow_damping_width
730       ENDIF
731    ENDIF
732
733!
734!-- Listing of 1D-profiles
735    WRITE ( io, 325 )  dt_dopr_listing
736    IF ( averaging_interval_pr /= 0.0 )  THEN
737       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
738    ENDIF
739
740!
741!-- DATA output
742    WRITE ( io, 330 )
743    IF ( averaging_interval_pr /= 0.0 )  THEN
744       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
745    ENDIF
746
747!
748!-- 1D-profiles
749    dopr_chr = 'Profile:'
750    IF ( dopr_n /= 0 )  THEN
751       WRITE ( io, 331 )
752
753       output_format = ''
754       IF ( netcdf_output )  THEN
755          IF ( netcdf_data_format == 1 )  THEN
756             output_format = 'NetCDF classic'
757          ELSE
758             output_format = 'NetCDF 64bit offset'
759          ENDIF
760       ENDIF
761       WRITE ( io, 344 )  output_format
762
763       DO  i = 1, dopr_n
764          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
765          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
766             WRITE ( io, 332 )  dopr_chr
767             dopr_chr = '       :'
768          ENDIF
769       ENDDO
770
771       IF ( dopr_chr /= '' )  THEN
772          WRITE ( io, 332 )  dopr_chr
773       ENDIF
774       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
775       IF ( skip_time_dopr /= 0.0 )  WRITE ( io, 339 )  skip_time_dopr
776    ENDIF
777
778!
779!-- 2D-arrays
780    DO  av = 0, 1
781
782       i = 1
783       do2d_xy = ''
784       do2d_xz = ''
785       do2d_yz = ''
786       DO  WHILE ( do2d(av,i) /= ' ' )
787
788          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
789          do2d_mode = do2d(av,i)(l-1:l)
790
791          SELECT CASE ( do2d_mode )
792             CASE ( 'xy' )
793                ll = LEN_TRIM( do2d_xy )
794                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
795             CASE ( 'xz' )
796                ll = LEN_TRIM( do2d_xz )
797                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
798             CASE ( 'yz' )
799                ll = LEN_TRIM( do2d_yz )
800                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
801          END SELECT
802
803          i = i + 1
804
805       ENDDO
806
807       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
808              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
809              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) )  .AND. &
810            ( netcdf_output  .OR.  iso2d_output ) )  THEN
811
812          IF (  av == 0 )  THEN
813             WRITE ( io, 334 )  ''
814          ELSE
815             WRITE ( io, 334 )  '(time-averaged)'
816          ENDIF
817
818          IF ( do2d_at_begin )  THEN
819             begin_chr = 'and at the start'
820          ELSE
821             begin_chr = ''
822          ENDIF
823
824          output_format = ''
825          IF ( netcdf_output )  THEN
826             IF ( netcdf_data_format == 1 )  THEN
827                output_format = 'NetCDF classic'
828             ELSEIF ( netcdf_data_format == 2 )  THEN
829                output_format = 'NetCDF 64bit offset'
830             ELSEIF ( netcdf_data_format == 3 )  THEN
831                output_format = 'NetCDF4/HDF5'
832             ELSEIF ( netcdf_data_format == 4 )  THEN
833                output_format = 'NetCDF4/HDF5 clasic'
834             ENDIF
835          ENDIF
836          IF ( iso2d_output )  THEN
837             IF ( netcdf_output )  THEN
838                output_format = TRIM( output_format ) // ' and iso2d'
839             ELSE
840                output_format = 'iso2d'
841             ENDIF
842          ENDIF
843          WRITE ( io, 344 )  output_format
844
845          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
846             i = 1
847             slices = '/'
848             coordinates = '/'
849!
850!--          Building strings with index and coordinate informations of the
851!--          slices
852             DO  WHILE ( section(i,1) /= -9999 )
853
854                WRITE (section_chr,'(I5)')  section(i,1)
855                section_chr = ADJUSTL( section_chr )
856                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
857
858                IF ( section(i,1) == -1 )  THEN
859                   WRITE (coor_chr,'(F10.1)')  -1.0
860                ELSE
861                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
862                ENDIF
863                coor_chr = ADJUSTL( coor_chr )
864                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
865
866                i = i + 1
867             ENDDO
868             IF ( av == 0 )  THEN
869                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
870                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
871                                   TRIM( coordinates )
872                IF ( skip_time_do2d_xy /= 0.0 )  THEN
873                   WRITE ( io, 339 )  skip_time_do2d_xy
874                ENDIF
875             ELSE
876                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
877                                   TRIM( begin_chr ), averaging_interval, &
878                                   dt_averaging_input, 'k', TRIM( slices ), &
879                                   TRIM( coordinates )
880                IF ( skip_time_data_output_av /= 0.0 )  THEN
881                   WRITE ( io, 339 )  skip_time_data_output_av
882                ENDIF
883             ENDIF
884
885          ENDIF
886
887          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
888             i = 1
889             slices = '/'
890             coordinates = '/'
891!
892!--          Building strings with index and coordinate informations of the
893!--          slices
894             DO  WHILE ( section(i,2) /= -9999 )
895
896                WRITE (section_chr,'(I5)')  section(i,2)
897                section_chr = ADJUSTL( section_chr )
898                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
899
900                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
901                coor_chr = ADJUSTL( coor_chr )
902                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
903
904                i = i + 1
905             ENDDO
906             IF ( av == 0 )  THEN
907                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
908                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
909                                   TRIM( coordinates )
910                IF ( skip_time_do2d_xz /= 0.0 )  THEN
911                   WRITE ( io, 339 )  skip_time_do2d_xz
912                ENDIF
913             ELSE
914                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
915                                   TRIM( begin_chr ), averaging_interval, &
916                                   dt_averaging_input, 'j', TRIM( slices ), &
917                                   TRIM( coordinates )
918                IF ( skip_time_data_output_av /= 0.0 )  THEN
919                   WRITE ( io, 339 )  skip_time_data_output_av
920                ENDIF
921             ENDIF
922          ENDIF
923
924          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
925             i = 1
926             slices = '/'
927             coordinates = '/'
928!
929!--          Building strings with index and coordinate informations of the
930!--          slices
931             DO  WHILE ( section(i,3) /= -9999 )
932
933                WRITE (section_chr,'(I5)')  section(i,3)
934                section_chr = ADJUSTL( section_chr )
935                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
936
937                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
938                coor_chr = ADJUSTL( coor_chr )
939                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
940
941                i = i + 1
942             ENDDO
943             IF ( av == 0 )  THEN
944                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
945                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
946                                   TRIM( coordinates )
947                IF ( skip_time_do2d_yz /= 0.0 )  THEN
948                   WRITE ( io, 339 )  skip_time_do2d_yz
949                ENDIF
950             ELSE
951                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
952                                   TRIM( begin_chr ), averaging_interval, &
953                                   dt_averaging_input, 'i', TRIM( slices ), &
954                                   TRIM( coordinates )
955                IF ( skip_time_data_output_av /= 0.0 )  THEN
956                   WRITE ( io, 339 )  skip_time_data_output_av
957                ENDIF
958             ENDIF
959          ENDIF
960
961       ENDIF
962
963    ENDDO
964
965!
966!-- 3d-arrays
967    DO  av = 0, 1
968
969       i = 1
970       do3d_chr = ''
971       DO  WHILE ( do3d(av,i) /= ' ' )
972
973          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
974          i = i + 1
975
976       ENDDO
977
978       IF ( do3d_chr /= '' )  THEN
979          IF ( av == 0 )  THEN
980             WRITE ( io, 336 )  ''
981          ELSE
982             WRITE ( io, 336 )  '(time-averaged)'
983          ENDIF
984
985          output_format = ''
986          IF ( netcdf_output )  THEN
987             IF ( netcdf_data_format == 1 )  THEN
988                output_format = 'NetCDF classic'
989             ELSEIF ( netcdf_data_format == 2 )  THEN
990                output_format = 'NetCDF 64bit offset'
991             ELSEIF ( netcdf_data_format == 3 )  THEN
992                output_format = 'NetCDF4/HDF5'
993             ELSEIF ( netcdf_data_format == 4 )  THEN
994                output_format = 'NetCDF4/HDF5 clasic'
995             ENDIF
996          ENDIF
997          IF ( avs_output )  THEN
998             IF ( netcdf_output )  THEN
999                output_format = TRIM( output_format ) // ' and avs'
1000             ELSE
1001                output_format = 'avs'
1002             ENDIF
1003          ENDIF
1004          WRITE ( io, 344 )  output_format
1005
1006          IF ( do3d_at_begin )  THEN
1007             begin_chr = 'and at the start'
1008          ELSE
1009             begin_chr = ''
1010          ENDIF
1011          IF ( av == 0 )  THEN
1012             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
1013                                zu(nz_do3d), nz_do3d
1014          ELSE
1015             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
1016                                TRIM( begin_chr ), averaging_interval, &
1017                                dt_averaging_input, zu(nz_do3d), nz_do3d
1018          ENDIF
1019
1020          IF ( do3d_compress )  THEN
1021             do3d_chr = ''
1022             i = 1
1023             DO WHILE ( do3d(av,i) /= ' ' )
1024
1025                SELECT CASE ( do3d(av,i) )
1026                   CASE ( 'u' )
1027                      j = 1
1028                   CASE ( 'v' )
1029                      j = 2
1030                   CASE ( 'w' )
1031                      j = 3
1032                   CASE ( 'p' )
1033                      j = 4
1034                   CASE ( 'pt' )
1035                      j = 5
1036                END SELECT
1037                WRITE ( prec, '(I1)' )  plot_3d_precision(j)%precision
1038                do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // &
1039                           ':' // prec // ','
1040                i = i + 1
1041
1042             ENDDO
1043             WRITE ( io, 338 )  do3d_chr
1044
1045          ENDIF
1046
1047          IF ( av == 0 )  THEN
1048             IF ( skip_time_do3d /= 0.0 )  THEN
1049                WRITE ( io, 339 )  skip_time_do3d
1050             ENDIF
1051          ELSE
1052             IF ( skip_time_data_output_av /= 0.0 )  THEN
1053                WRITE ( io, 339 )  skip_time_data_output_av
1054             ENDIF
1055          ENDIF
1056
1057       ENDIF
1058
1059    ENDDO
1060
1061!
1062!-- masked arrays
1063    IF ( masks > 0 )  WRITE ( io, 345 )  &
1064         mask_scale_x, mask_scale_y, mask_scale_z
1065    DO  mid = 1, masks
1066       DO  av = 0, 1
1067
1068          i = 1
1069          domask_chr = ''
1070          DO  WHILE ( domask(mid,av,i) /= ' ' )
1071             domask_chr = TRIM( domask_chr ) // ' ' //  &
1072                          TRIM( domask(mid,av,i) ) // ','
1073             i = i + 1
1074          ENDDO
1075
1076          IF ( domask_chr /= '' )  THEN
1077             IF ( av == 0 )  THEN
1078                WRITE ( io, 346 )  '', mid
1079             ELSE
1080                WRITE ( io, 346 )  ' (time-averaged)', mid
1081             ENDIF
1082
1083             output_format = ''
1084             IF ( netcdf_output )  THEN
1085                IF ( netcdf_data_format == 1 )  THEN
1086                   output_format = 'NetCDF classic'
1087                ELSEIF ( netcdf_data_format == 2 )  THEN
1088                   output_format = 'NetCDF 64bit offset'
1089                ELSEIF ( netcdf_data_format == 3 )  THEN
1090                   output_format = 'NetCDF4/HDF5'
1091                ELSEIF ( netcdf_data_format == 4 )  THEN
1092                   output_format = 'NetCDF4/HDF5 clasic'
1093                ENDIF
1094             ENDIF
1095             WRITE ( io, 344 )  output_format
1096
1097             IF ( av == 0 )  THEN
1098                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1099             ELSE
1100                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1101                                   averaging_interval, dt_averaging_input
1102             ENDIF
1103
1104             IF ( av == 0 )  THEN
1105                IF ( skip_time_domask(mid) /= 0.0 )  THEN
1106                   WRITE ( io, 339 )  skip_time_domask(mid)
1107                ENDIF
1108             ELSE
1109                IF ( skip_time_data_output_av /= 0.0 )  THEN
1110                   WRITE ( io, 339 )  skip_time_data_output_av
1111                ENDIF
1112             ENDIF
1113!
1114!--          output locations
1115             DO  dim = 1, 3
1116                IF ( mask(mid,dim,1) >= 0.0 )  THEN
1117                   count = 0
1118                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0 )
1119                      count = count + 1
1120                   ENDDO
1121                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1122                                      mask(mid,dim,:count)
1123                ELSEIF ( mask_loop(mid,dim,1) < 0.0 .AND.  &
1124                         mask_loop(mid,dim,2) < 0.0 .AND.  &
1125                         mask_loop(mid,dim,3) == 0.0 )  THEN
1126                   WRITE ( io, 350 )  dir(dim), dir(dim)
1127                ELSEIF ( mask_loop(mid,dim,3) == 0.0 )  THEN
1128                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1129                                      mask_loop(mid,dim,1:2)
1130                ELSE
1131                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1132                                      mask_loop(mid,dim,1:3)
1133                ENDIF
1134             ENDDO
1135          ENDIF
1136
1137       ENDDO
1138    ENDDO
1139
1140!
1141!-- Timeseries
1142    IF ( dt_dots /= 9999999.9 )  THEN
1143       WRITE ( io, 340 )
1144
1145       output_format = ''
1146       IF ( netcdf_output )  THEN
1147          IF ( netcdf_data_format == 1 )  THEN
1148             output_format = 'NetCDF classic'
1149          ELSE
1150             output_format = 'NetCDF 64bit offset'
1151          ENDIF
1152       ENDIF
1153       WRITE ( io, 344 )  output_format
1154       WRITE ( io, 341 )  dt_dots
1155    ENDIF
1156
1157#if defined( __dvrp_graphics )
1158!
1159!-- Dvrp-output
1160    IF ( dt_dvrp /= 9999999.9 )  THEN
1161       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1162                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1163       i = 1
1164       l = 0
1165       m = 0
1166       DO WHILE ( mode_dvrp(i) /= ' ' )
1167          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
1168             READ ( mode_dvrp(i), '(10X,I2)' )  j
1169             l = l + 1
1170             IF ( do3d(0,j) /= ' ' )  THEN
1171                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1172                                   isosurface_color(:,l)
1173             ENDIF
1174          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
1175             READ ( mode_dvrp(i), '(6X,I2)' )  j
1176             m = m + 1
1177             IF ( do2d(0,j) /= ' ' )  THEN
1178                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1179                                   slicer_range_limits_dvrp(:,m)
1180             ENDIF
1181          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
1182             WRITE ( io, 363 )  dvrp_psize
1183             IF ( particle_dvrpsize /= 'none' )  THEN
1184                WRITE ( io, 364 )  'size', TRIM( particle_dvrpsize ), &
1185                                   dvrpsize_interval
1186             ENDIF
1187             IF ( particle_color /= 'none' )  THEN
1188                WRITE ( io, 364 )  'color', TRIM( particle_color ), &
1189                                   color_interval
1190             ENDIF
1191          ENDIF
1192          i = i + 1
1193       ENDDO
1194
1195       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1196                          superelevation_y, superelevation, clip_dvrp_l, &
1197                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1198
1199       IF ( TRIM( topography ) /= 'flat' )  THEN
1200          WRITE ( io, 366 )  topography_color
1201          IF ( cluster_size > 1 )  THEN
1202             WRITE ( io, 367 )  cluster_size
1203          ENDIF
1204       ENDIF
1205
1206    ENDIF
1207#endif
1208
1209#if defined( __spectra )
1210!
1211!-- Spectra output
1212    IF ( dt_dosp /= 9999999.9 ) THEN
1213       WRITE ( io, 370 )
1214
1215       output_format = ''
1216       IF ( netcdf_output )  THEN
1217          IF ( netcdf_data_format == 1 )  THEN
1218             output_format = 'NetCDF classic'
1219          ELSE
1220             output_format = 'NetCDF 64bit offset'
1221          ENDIF
1222       ENDIF
1223       WRITE ( io, 344 )  output_format
1224       WRITE ( io, 371 )  dt_dosp
1225       IF ( skip_time_dosp /= 0.0 )  WRITE ( io, 339 )  skip_time_dosp
1226       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
1227                          ( spectra_direction(i), i = 1,10 ),  &
1228                          ( comp_spectra_level(i), i = 1,100 ), &
1229                          ( plot_spectra_level(i), i = 1,100 ), &
1230                          averaging_interval_sp, dt_averaging_input_pr
1231    ENDIF
1232#endif
1233
1234    WRITE ( io, 99 )
1235
1236!
1237!-- Physical quantities
1238    WRITE ( io, 400 )
1239
1240!
1241!-- Geostrophic parameters
1242    WRITE ( io, 410 )  omega, phi, f, fs
1243
1244!
1245!-- Other quantities
1246    WRITE ( io, 411 )  g
1247    IF ( use_reference )  THEN
1248       IF ( ocean )  THEN
1249          WRITE ( io, 412 )  prho_reference
1250       ELSE
1251          WRITE ( io, 413 )  pt_reference
1252       ENDIF
1253    ENDIF
1254
1255!
1256!-- Cloud physics parameters
1257    IF ( cloud_physics ) THEN
1258       WRITE ( io, 415 )
1259       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
1260    ENDIF
1261
1262!-- Profile of the geostrophic wind (component ug)
1263!-- Building output strings
1264    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
1265    gradients = '------'
1266    slices = '     0'
1267    coordinates = '   0.0'
1268    i = 1
1269    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
1270     
1271       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
1272       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
1273
1274       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
1275       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1276
1277       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
1278       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1279
1280       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
1281       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1282
1283       IF ( i == 10 )  THEN
1284          EXIT
1285       ELSE
1286          i = i + 1
1287       ENDIF
1288
1289    ENDDO
1290
1291    WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
1292                       TRIM( gradients ), TRIM( slices )
1293
1294!-- Profile of the geostrophic wind (component vg)
1295!-- Building output strings
1296    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
1297    gradients = '------'
1298    slices = '     0'
1299    coordinates = '   0.0'
1300    i = 1
1301    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
1302
1303       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
1304       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
1305
1306       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
1307       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1308
1309       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
1310       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1311
1312       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
1313       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1314
1315       IF ( i == 10 )  THEN
1316          EXIT
1317       ELSE
1318          i = i + 1
1319       ENDIF
1320 
1321    ENDDO
1322
1323    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
1324                       TRIM( gradients ), TRIM( slices )
1325
1326!
1327!-- Initial wind profiles
1328    IF ( u_profile(1) /= 9999999.9 )  WRITE ( io, 427 )
1329
1330!
1331!-- Initial temperature profile
1332!-- Building output strings, starting with surface temperature
1333    WRITE ( temperatures, '(F6.2)' )  pt_surface
1334    gradients = '------'
1335    slices = '     0'
1336    coordinates = '   0.0'
1337    i = 1
1338    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1339
1340       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1341       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1342
1343       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1344       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1345
1346       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1347       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1348
1349       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1350       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1351
1352       IF ( i == 10 )  THEN
1353          EXIT
1354       ELSE
1355          i = i + 1
1356       ENDIF
1357
1358    ENDDO
1359
1360    WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1361                       TRIM( gradients ), TRIM( slices )
1362
1363!
1364!-- Initial humidity profile
1365!-- Building output strings, starting with surface humidity
1366    IF ( humidity  .OR.  passive_scalar )  THEN
1367       WRITE ( temperatures, '(E8.1)' )  q_surface
1368       gradients = '--------'
1369       slices = '       0'
1370       coordinates = '     0.0'
1371       i = 1
1372       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1373         
1374          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1375          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1376
1377          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1378          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1379         
1380          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1381          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1382         
1383          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1384          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1385
1386          IF ( i == 10 )  THEN
1387             EXIT
1388          ELSE
1389             i = i + 1
1390          ENDIF
1391
1392       ENDDO
1393
1394       IF ( humidity )  THEN
1395          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1396                             TRIM( gradients ), TRIM( slices )
1397       ELSE
1398          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1399                             TRIM( gradients ), TRIM( slices )
1400       ENDIF
1401    ENDIF
1402
1403!
1404!-- Initial salinity profile
1405!-- Building output strings, starting with surface salinity
1406    IF ( ocean )  THEN
1407       WRITE ( temperatures, '(F6.2)' )  sa_surface
1408       gradients = '------'
1409       slices = '     0'
1410       coordinates = '   0.0'
1411       i = 1
1412       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1413
1414          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1415          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1416
1417          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1418          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1419
1420          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1421          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1422
1423          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1424          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1425
1426          IF ( i == 10 )  THEN
1427             EXIT
1428          ELSE
1429             i = i + 1
1430          ENDIF
1431
1432       ENDDO
1433
1434       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1435                          TRIM( gradients ), TRIM( slices )
1436    ENDIF
1437
1438!
1439!-- Profile for the large scale vertial velocity
1440!-- Building output strings, starting with surface value
1441    IF ( large_scale_subsidence )  THEN
1442       temperatures = '   0.0'
1443       gradients = '------'
1444       slices = '     0'
1445       coordinates = '   0.0'
1446       i = 1
1447       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
1448
1449          WRITE (coor_chr,'(E10.2,7X)')  &
1450                                w_subs(subs_vertical_gradient_level_i(i))
1451          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1452
1453          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
1454          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1455
1456          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
1457          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1458
1459          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
1460          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1461
1462          IF ( i == 10 )  THEN
1463             EXIT
1464          ELSE
1465             i = i + 1
1466          ENDIF
1467
1468       ENDDO
1469
1470       WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
1471                          TRIM( gradients ), TRIM( slices )
1472    ENDIF
1473
1474!
1475!-- Cloud physcis parameters / quantities / numerical methods
1476    WRITE ( io, 430 )
1477    IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets)  THEN
1478       WRITE ( io, 431 )
1479    ELSEIF ( humidity  .AND.  cloud_physics )  THEN
1480       WRITE ( io, 432 )
1481       IF ( radiation )      WRITE ( io, 132 )
1482       IF ( precipitation )  WRITE ( io, 133 )
1483    ELSEIF ( humidity  .AND.  cloud_droplets )  THEN
1484       WRITE ( io, 433 )
1485       IF ( curvature_solution_effects )  WRITE ( io, 434 )
1486       IF ( collision_kernel /= 'none' )  THEN
1487          WRITE ( io, 435 )  TRIM( collision_kernel )
1488          IF ( collision_kernel(6:9) == 'fast' )  THEN
1489             WRITE ( io, 436 )  radius_classes, dissipation_classes
1490          ENDIF
1491       ELSE
1492          WRITE ( io, 437 )
1493       ENDIF
1494    ENDIF
1495
1496!
1497!-- LES / turbulence parameters
1498    WRITE ( io, 450 )
1499
1500!--
1501! ... LES-constants used must still be added here
1502!--
1503    IF ( constant_diffusion )  THEN
1504       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1505                          prandtl_number
1506    ENDIF
1507    IF ( .NOT. constant_diffusion)  THEN
1508       IF ( e_init > 0.0 )  WRITE ( io, 455 )  e_init
1509       IF ( e_min > 0.0 )  WRITE ( io, 454 )  e_min
1510       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1511       IF ( adjust_mixing_length  .AND.  prandtl_layer )  WRITE ( io, 452 )
1512    ENDIF
1513
1514!
1515!-- Special actions during the run
1516    WRITE ( io, 470 )
1517    IF ( create_disturbances )  THEN
1518       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1519                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1520                          zu(disturbance_level_ind_t), disturbance_level_ind_t
1521       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
1522          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1523       ELSE
1524          WRITE ( io, 473 )  disturbance_energy_limit
1525       ENDIF
1526       WRITE ( io, 474 )  TRIM( random_generator )
1527    ENDIF
1528    IF ( pt_surface_initial_change /= 0.0 )  THEN
1529       WRITE ( io, 475 )  pt_surface_initial_change
1530    ENDIF
1531    IF ( humidity  .AND.  q_surface_initial_change /= 0.0 )  THEN
1532       WRITE ( io, 476 )  q_surface_initial_change       
1533    ENDIF
1534    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0 )  THEN
1535       WRITE ( io, 477 )  q_surface_initial_change       
1536    ENDIF
1537
1538    IF ( particle_advection )  THEN
1539!
1540!--    Particle attributes
1541       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1542                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
1543                          end_time_prel, dt_sort_particles
1544       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1545       IF ( random_start_position )  WRITE ( io, 481 )
1546       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1547       WRITE ( io, 495 )  total_number_of_particles
1548       IF ( use_particle_tails  .AND.  maximum_number_of_tailpoints /= 0 )  THEN
1549          WRITE ( io, 483 )  maximum_number_of_tailpoints
1550          IF ( minimum_tailpoint_distance /= 0 )  THEN
1551             WRITE ( io, 484 )  total_number_of_tails,      &
1552                                minimum_tailpoint_distance, &
1553                                maximum_tailpoint_age
1554          ENDIF
1555       ENDIF
1556       IF ( dt_write_particle_data /= 9999999.9 )  THEN
1557          WRITE ( io, 485 )  dt_write_particle_data
1558          output_format = ''
1559          IF ( netcdf_output )  THEN
1560             IF ( netcdf_data_format > 1 )  THEN
1561                output_format = 'netcdf (64 bit offset) and binary'
1562             ELSE
1563                output_format = 'netcdf and binary'
1564             ENDIF
1565          ELSE
1566             output_format = 'binary'
1567          ENDIF
1568          WRITE ( io, 344 )  output_format
1569       ENDIF
1570       IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1571       IF ( write_particle_statistics )  WRITE ( io, 486 )
1572
1573       WRITE ( io, 487 )  number_of_particle_groups
1574
1575       DO  i = 1, number_of_particle_groups
1576          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9 )  THEN
1577             WRITE ( io, 490 )  i, 0.0
1578             WRITE ( io, 492 )
1579          ELSE
1580             WRITE ( io, 490 )  i, radius(i)
1581             IF ( density_ratio(i) /= 0.0 )  THEN
1582                WRITE ( io, 491 )  density_ratio(i)
1583             ELSE
1584                WRITE ( io, 492 )
1585             ENDIF
1586          ENDIF
1587          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1588                             pdx(i), pdy(i), pdz(i)
1589          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
1590       ENDDO
1591
1592    ENDIF
1593
1594
1595!
1596!-- Parameters of 1D-model
1597    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1598       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1599                          mixing_length_1d, dissipation_1d
1600       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1601          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1602       ENDIF
1603    ENDIF
1604
1605!
1606!-- User-defined informations
1607    CALL user_header( io )
1608
1609    WRITE ( io, 99 )
1610
1611!
1612!-- Write buffer contents to disc immediately
1613    CALL local_flush( io )
1614
1615!
1616!-- Here the FORMATs start
1617
1618 99 FORMAT (1X,78('-'))
1619100 FORMAT (/1X,'***************************',9X,42('-')/        &
1620            1X,'* ',A,' *',9X,A/                               &
1621            1X,'***************************',9X,42('-'))
1622101 FORMAT (37X,'coupled run using MPI-',I1,': ',A/ &
1623            37X,42('-'))
1624102 FORMAT (/' Date:              ',A8,9X,'Run:       ',A20/      &
1625            ' Time:              ',A8,9X,'Run-No.:   ',I2.2/     &
1626            ' Run on host:     ',A10)
1627#if defined( __parallel )
1628103 FORMAT (' Number of PEs:',8X,I5,9X,'Processor grid (x,y): (',I3,',',I3, &
1629              ')',1X,A)
1630104 FORMAT (' Number of PEs:',8X,I5,9X,'Tasks:',I4,'   threads per task:',I4/ &
1631              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
1632105 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1633106 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
1634            37X,'because the job is running on an SMP-cluster')
1635107 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
1636108 FORMAT (37X,'Max. # of parallel I/O streams is ',I5)
1637#endif
1638110 FORMAT (/' Numerical Schemes:'/ &
1639             ' -----------------'/)
1640111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1641112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1642            '     Iterations (initial/other): ',I3,'/',I3,'  omega = ',F5.3)
1643113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1644                  ' or Upstream')
1645114 FORMAT (' --> Momentum advection via Upstream-Spline-Scheme')
1646115 FORMAT ('     Tendencies are smoothed via Long-Filter with factor ',F5.3) 
1647116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1648                  ' or Upstream')
1649117 FORMAT (' --> Scalar advection via Upstream-Spline-Scheme')
1650118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1651119 FORMAT (' --> Galilei-Transform applied to horizontal advection', &
1652            '     Translation velocity = ',A/ &
1653            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1654120 FORMAT (' --> Time differencing scheme: leapfrog only (no euler in case', &
1655                  ' of timestep changes)')
1656121 FORMAT (' --> Time differencing scheme: leapfrog + euler in case of', &
1657                  ' timestep changes')
1658122 FORMAT (' --> Time differencing scheme: ',A)
1659123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
1660            '     maximum damping coefficient: ',F5.3, ' 1/s')
1661124 FORMAT ('     Spline-overshoots are being suppressed')
1662125 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1663                  ' of'/                                                       &
1664            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1665126 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1666                  ' of'/                                                       &
1667            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1668127 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1669            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1670128 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1671            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1672129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1673130 FORMAT (' --> Additional prognostic equation for the total water content')
1674131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
1675                  F6.2, ' K assumed')
1676132 FORMAT ('     Parameterization of long-wave radiation processes via'/ &
1677            '     effective emissivity scheme')
1678133 FORMAT ('     Precipitation parameterization via Kessler-Scheme')
1679134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1680135 FORMAT (' --> Solve perturbation pressure via multigrid method (', &
1681                  A,'-cycle)'/ &
1682            '     number of grid levels:                   ',I2/ &
1683            '     Gauss-Seidel red/black iterations:       ',I2)
1684136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1685                  I3,')')
1686137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1687            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1688                  I3,')'/ &
1689            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1690                  I3,')')
1691138 FORMAT ('     Using hybrid version for 1d-domain-decomposition')
1692139 FORMAT (' --> Loop optimization method: ',A)
1693140 FORMAT ('     maximum residual allowed:                ',E10.3)
1694141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1695142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1696                  'step')
1697143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1698                  'kinetic energy')
1699144 FORMAT ('     masking method is used')
1700150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
1701                  'conserved'/ &
1702            '     using the ',A,' mode')
1703151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
1704152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
1705           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
1706           /'     starting from dp_level_b =', F8.3, 'm', A /)
1707153 FORMAT (' --> Large-scale vertical motion is used in the ', &
1708                  'prognostic equation for')
1709154 FORMAT ('     the potential temperature')
1710200 FORMAT (//' Run time and time step information:'/ &
1711             ' ----------------------------------'/)
1712201 FORMAT ( ' Timestep:          variable     maximum value: ',F6.3,' s', &
1713             '    CFL-factor: ',F4.2)
1714202 FORMAT ( ' Timestep:       dt = ',F6.3,' s'/)
1715203 FORMAT ( ' Start time:       ',F9.3,' s'/ &
1716             ' End time:         ',F9.3,' s')
1717204 FORMAT ( A,F9.3,' s')
1718205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
1719206 FORMAT (/' Time reached:     ',F9.3,' s'/ &
1720             ' CPU-time used:    ',F9.3,' s     per timestep:               ', &
1721               '  ',F9.3,' s'/                                                 &
1722             '                                   per second of simulated tim', &
1723               'e: ',F9.3,' s')
1724207 FORMAT ( A/' Coupling start time:',F9.3,' s')
1725250 FORMAT (//' Computational grid and domain size:'/ &
1726              ' ----------------------------------'// &
1727              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1728              ' m    dz =    ',F7.3,' m'/ &
1729              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1730              ' m  z(u) = ',F10.3,' m'/)
1731252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
1732              ' factor: ',F5.3/ &
1733            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1734254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1735            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1736255 FORMAT (' Subdomains have equal size')
1737256 FORMAT (' Subdomains at the upper edges of the virtual processor grid ', &
1738              'have smaller sizes'/                                          &
1739            ' Size of smallest subdomain:    (  ',I4,',   ',I4,',   ',I4,')')
1740260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1741             ' degrees')
1742270 FORMAT (//' Topography informations:'/ &
1743              ' -----------------------'// &
1744              1X,'Topography: ',A)
1745271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1746              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1747                ' / ',I4)
1748272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
1749              ' direction' / &
1750              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
1751              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
1752278 FORMAT (' Topography grid definition convention:'/ &
1753            ' cell edge (staggered grid points'/  &
1754            ' (u in x-direction, v in y-direction))' /)
1755279 FORMAT (' Topography grid definition convention:'/ &
1756            ' cell center (scalar grid points)' /)
1757280 FORMAT (//' Vegetation canopy (drag) model:'/ &
1758              ' ------------------------------'// &
1759              ' Canopy mode: ', A / &
1760              ' Canopy top: ',I4 / &
1761              ' Leaf drag coefficient: ',F6.2 /)
1762281 FORMAT (/ ' Scalar_exchange_coefficient: ',F6.2 / &
1763              ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /)
1764282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s')
1765283 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
1766              ' Height:              ',A,'  m'/ &
1767              ' Leaf area density:   ',A,'  m**2/m**3'/ &
1768              ' Gradient:            ',A,'  m**2/m**4'/ &
1769              ' Gridpoint:           ',A)
1770               
1771300 FORMAT (//' Boundary conditions:'/ &
1772             ' -------------------'// &
1773             '                     p                    uv             ', &
1774             '                   pt'// &
1775             ' B. bound.: ',A/ &
1776             ' T. bound.: ',A)
1777301 FORMAT (/'                     ',A// &
1778             ' B. bound.: ',A/ &
1779             ' T. bound.: ',A)
1780303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1781304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1782305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
1783               'computational u,v-level:'// &
1784             '       zp = ',F6.2,' m   z0 = ',F6.4,' m   z0h = ',F7.5,&
1785             ' m   kappa = ',F4.2/ &
1786             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1787306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
1788307 FORMAT ('       Heatflux has a random normal distribution')
1789308 FORMAT ('       Predefined surface temperature')
1790309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
1791310 FORMAT (//'    1D-Model:'// &
1792             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1793311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1794312 FORMAT ('       Predefined surface humidity')
1795313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1796314 FORMAT ('       Predefined scalar value at the surface')
1797315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
1798316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1799                    'atmosphere model')
1800317 FORMAT (//' Lateral boundaries:'/ &
1801            '       left/right:  ',A/    &
1802            '       north/south: ',A)
1803318 FORMAT (/'       pt damping layer width = ',F7.2,' m, pt ', &
1804                    'damping factor = ',F6.4)
1805319 FORMAT ('       turbulence recycling at inflow switched on'/ &
1806            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1807            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
1808320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
1809            '                                          v: ',F9.6,' m**2/s**2')
1810325 FORMAT (//' List output:'/ &
1811             ' -----------'//  &
1812            '    1D-Profiles:'/    &
1813            '       Output every             ',F8.2,' s')
1814326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
1815            '       Averaging input every    ',F8.2,' s')
1816330 FORMAT (//' Data output:'/ &
1817             ' -----------'/)
1818331 FORMAT (/'    1D-Profiles:')
1819332 FORMAT (/'       ',A)
1820333 FORMAT ('       Output every             ',F8.2,' s',/ &
1821            '       Time averaged over       ',F8.2,' s'/ &
1822            '       Averaging input every    ',F8.2,' s')
1823334 FORMAT (/'    2D-Arrays',A,':')
1824335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1825            '       Output every             ',F8.2,' s  ',A/ &
1826            '       Cross sections at ',A1,' = ',A/ &
1827            '       scalar-coordinates:   ',A,' m'/)
1828336 FORMAT (/'    3D-Arrays',A,':')
1829337 FORMAT (/'       Arrays: ',A/ &
1830            '       Output every             ',F8.2,' s  ',A/ &
1831            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1832338 FORMAT ('       Compressed data output'/ &
1833            '       Decimal precision: ',A/)
1834339 FORMAT ('       No output during initial ',F8.2,' s')
1835340 FORMAT (/'    Time series:')
1836341 FORMAT ('       Output every             ',F8.2,' s'/)
1837342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1838            '       Output every             ',F8.2,' s  ',A/ &
1839            '       Time averaged over       ',F8.2,' s'/ &
1840            '       Averaging input every    ',F8.2,' s'/ &
1841            '       Cross sections at ',A1,' = ',A/ &
1842            '       scalar-coordinates:   ',A,' m'/)
1843343 FORMAT (/'       Arrays: ',A/ &
1844            '       Output every             ',F8.2,' s  ',A/ &
1845            '       Time averaged over       ',F8.2,' s'/ &
1846            '       Averaging input every    ',F8.2,' s'/ &
1847            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1848344 FORMAT ('       Output format: ',A/)
1849345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
1850            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
1851            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
1852            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
1853346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
1854347 FORMAT ('       Variables: ',A/ &
1855            '       Output every             ',F8.2,' s')
1856348 FORMAT ('       Variables: ',A/ &
1857            '       Output every             ',F8.2,' s'/ &
1858            '       Time averaged over       ',F8.2,' s'/ &
1859            '       Averaging input every    ',F8.2,' s')
1860349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1861            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
1862            13('       ',8(F8.2,',')/) )
1863350 FORMAT (/'       Output locations in ',A,'-direction: ', &
1864            'all gridpoints along ',A,'-direction (default).' )
1865351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1866            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
1867            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
1868#if defined( __dvrp_graphics )
1869360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
1870            '       Output every      ',F7.1,' s'/ &
1871            '       Output mode:      ',A/ &
1872            '       Host / User:      ',A,' / ',A/ &
1873            '       Directory:        ',A// &
1874            '       The sequence contains:')
1875361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
1876            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
1877362 FORMAT (/'       Slicer plane ',A/ &
1878            '       Slicer limits: [',F6.2,',',F6.2,']')
1879363 FORMAT (/'       Particles'/ &
1880            '          particle size:  ',F7.2,' m')
1881364 FORMAT ('          particle ',A,' controlled by "',A,'" with interval [', &
1882                       F6.2,',',F6.2,']')
1883365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
1884            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
1885                     ')'/ &
1886            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
1887            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
1888366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
1889367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
1890#endif
1891#if defined( __spectra )
1892370 FORMAT ('    Spectra:')
1893371 FORMAT ('       Output every ',F7.1,' s'/)
1894372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
1895            '       Directions: ', 10(A5,',')/                         &
1896            '       height levels  k = ', 20(I3,',')/                  &
1897            '                          ', 20(I3,',')/                  &
1898            '                          ', 20(I3,',')/                  &
1899            '                          ', 20(I3,',')/                  &
1900            '                          ', 19(I3,','),I3,'.'/           &
1901            '       height levels selected for standard plot:'/        &
1902            '                      k = ', 20(I3,',')/                  &
1903            '                          ', 20(I3,',')/                  &
1904            '                          ', 20(I3,',')/                  &
1905            '                          ', 20(I3,',')/                  &
1906            '                          ', 19(I3,','),I3,'.'/           &
1907            '       Time averaged over ', F7.1, ' s,' /                &
1908            '       Profiles for the time averaging are taken every ', &
1909                    F6.1,' s')
1910#endif
1911400 FORMAT (//' Physical quantities:'/ &
1912              ' -------------------'/)
1913410 FORMAT ('    Angular velocity    :   omega = ',E9.3,' rad/s'/  &
1914            '    Geograph. latitude  :   phi   = ',F4.1,' degr'/   &
1915            '    Coriolis parameter  :   f     = ',F9.6,' 1/s'/    &
1916            '                            f*    = ',F9.6,' 1/s')
1917411 FORMAT (/'    Gravity             :   g     = ',F4.1,' m/s**2')
1918412 FORMAT (/'    Reference density in buoyancy terms: ',F8.3,' kg/m**3')
1919413 FORMAT (/'    Reference temperature in buoyancy terms: ',F8.4,' K')
1920415 FORMAT (/'    Cloud physics parameters:'/ &
1921             '    ------------------------'/)
1922416 FORMAT ('        Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
1923            '        Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
1924            '        Density of air     :   rho_0 = ',F5.3,' kg/m**3'/  &
1925            '        Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
1926            '        Vapourization heat :   L_v   = ',E8.2,' J/kg')
1927420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1928            '       Height:        ',A,'  m'/ &
1929            '       Temperature:   ',A,'  K'/ &
1930            '       Gradient:      ',A,'  K/100m'/ &
1931            '       Gridpoint:     ',A)
1932421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1933            '       Height:      ',A,'  m'/ &
1934            '       Humidity:    ',A,'  kg/kg'/ &
1935            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1936            '       Gridpoint:   ',A)
1937422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1938            '       Height:                  ',A,'  m'/ &
1939            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1940            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1941            '       Gridpoint:               ',A)
1942423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1943            '       Height:      ',A,'  m'/ &
1944            '       ug:          ',A,'  m/s'/ &
1945            '       Gradient:    ',A,'  1/100s'/ &
1946            '       Gridpoint:   ',A)
1947424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1948            '       Height:      ',A,'  m'/ &
1949            '       vg:          ',A,'  m/s'/ &
1950            '       Gradient:    ',A,'  1/100s'/ &
1951            '       Gridpoint:   ',A)
1952425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
1953            '       Height:     ',A,'  m'/ &
1954            '       Salinity:   ',A,'  psu'/ &
1955            '       Gradient:   ',A,'  psu/100m'/ &
1956            '       Gridpoint:  ',A)
1957426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
1958            '       Height:      ',A,'  m'/ &
1959            '       w_subs:      ',A,'  m/s'/ &
1960            '       Gradient:    ',A,'  (m/s)/100m'/ &
1961            '       Gridpoint:   ',A)
1962427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
1963                  ' profiles')
1964430 FORMAT (//' Cloud physics quantities / methods:'/ &
1965              ' ----------------------------------'/)
1966431 FORMAT ('    Humidity is treated as purely passive scalar (no condensati', &
1967                 'on)')
1968432 FORMAT ('    Bulk scheme with liquid water potential temperature and'/ &
1969            '    total water content is used.'/ &
1970            '    Condensation is parameterized via 0% - or 100% scheme.')
1971433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
1972                 'icle model')
1973434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
1974                 ' droplets < 1.0E-6 m')
1975435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
1976436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
1977                    'are used'/ &
1978            '          number of radius classes:       ',I3,'    interval ', &
1979                       '[1.0E-6,2.0E-4] m'/ &
1980            '          number of dissipation classes:   ',I2,'    interval ', &
1981                       '[0,1000] cm**2/s**3')
1982437 FORMAT ('    Droplet collision is switched off')
1983450 FORMAT (//' LES / Turbulence quantities:'/ &
1984              ' ---------------------------'/)
1985451 FORMAT ('    Diffusion coefficients are constant:'/ &
1986            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1987452 FORMAT ('    Mixing length is limited to the Prandtl mixing lenth.')
1988453 FORMAT ('    Mixing length is limited to ',F4.2,' * z')
1989454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
1990455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
1991470 FORMAT (//' Actions during the simulation:'/ &
1992              ' -----------------------------'/)
1993471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1994            '    Disturbance amplitude           :     ',F4.2, ' m/s'/       &
1995            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1996            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
1997472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1998                 ' to i/j =',I4)
1999473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
2000                 1X,F5.3, ' m**2/s**2')
2001474 FORMAT ('    Random number generator used    : ',A/)
2002475 FORMAT ('    The surface temperature is increased (or decreased, ', &
2003                 'respectively, if'/ &
2004            '    the value is negative) by ',F5.2,' K at the beginning of the',&
2005                 ' 3D-simulation'/)
2006476 FORMAT ('    The surface humidity is increased (or decreased, ',&
2007                 'respectively, if the'/ &
2008            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
2009                 ' the 3D-simulation'/)
2010477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
2011                 'respectively, if the'/ &
2012            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
2013                 ' the 3D-simulation'/)
2014480 FORMAT ('    Particles:'/ &
2015            '    ---------'// &
2016            '       Particle advection is active (switched on at t = ', F7.1, &
2017                    ' s)'/ &
2018            '       Start of new particle generations every  ',F6.1,' s'/ &
2019            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
2020            '                            bottom:     ', A, ' top:         ', A/&
2021            '       Maximum particle age:                 ',F9.1,' s'/ &
2022            '       Advection stopped at t = ',F9.1,' s'/ &
2023            '       Particles are sorted every ',F9.1,' s'/)
2024481 FORMAT ('       Particles have random start positions'/)
2025482 FORMAT ('          Particles are advected only horizontally'/)
2026483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
2027484 FORMAT ('            Number of tails of the total domain: ',I10/ &
2028            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
2029            '            Maximum age of the end of the tail:  ',F8.2,' s')
2030485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
2031486 FORMAT ('       Particle statistics are written on file'/)
2032487 FORMAT ('       Number of particle groups: ',I2/)
2033488 FORMAT ('       SGS velocity components are used for particle advection'/ &
2034            '          minimum timestep for advection: ', F7.5/)
2035489 FORMAT ('       Number of particles simultaneously released at each ', &
2036                    'point: ', I5/)
2037490 FORMAT ('       Particle group ',I2,':'/ &
2038            '          Particle radius: ',E10.3, 'm')
2039491 FORMAT ('          Particle inertia is activated'/ &
2040            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
2041492 FORMAT ('          Particles are advected only passively (no inertia)'/)
2042493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
2043            '                                         y:',F8.1,' - ',F8.1,' m'/&
2044            '                                         z:',F8.1,' - ',F8.1,' m'/&
2045            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
2046                       ' m  dz = ',F8.1,' m'/)
2047494 FORMAT ('       Output of particle time series in NetCDF format every ', &
2048                    F8.2,' s'/)
2049495 FORMAT ('       Number of particles in total domain: ',I10/)
2050500 FORMAT (//' 1D-Model parameters:'/                           &
2051              ' -------------------'//                           &
2052            '    Simulation time:                   ',F8.1,' s'/ &
2053            '    Run-controll output every:         ',F8.1,' s'/ &
2054            '    Vertical profile output every:     ',F8.1,' s'/ &
2055            '    Mixing length calculation:         ',A/         &
2056            '    Dissipation calculation:           ',A/)
2057502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
2058503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2059504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
2060
2061
2062 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.