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

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

last commit documented

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