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

Last change on this file since 927 was 927, checked in by raasch, 12 years ago

masking method can be switched on for mg-solver using inipar parameter masking_method

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