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

Last change on this file since 138 was 138, checked in by letzel, 16 years ago

Plant canopy model of Watanabe (2004,BLM 112,307-341) added.

  • Property svn:keywords set to Id
File size: 60.8 KB
Line 
1 SUBROUTINE header
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
7! Allow two instead of one digit to specify isosurface and slicer variables.
8! Output of sorting frequency of particles
9!
10! Former revisions:
11! -----------------
12! $Id: header.f90 138 2007-11-28 10:03:58Z letzel $
13!
14! 108 2007-08-24 15:10:38Z letzel
15! Output of informations for coupled model runs (boundary conditions etc.)
16! + output of momentumfluxes at the top boundary
17! Rayleigh damping for ocean, e_init
18!
19! 97 2007-06-21 08:23:15Z raasch
20! Adjustments for the ocean version.
21! use_pt_reference renamed use_reference
22!
23! 87 2007-05-22 15:46:47Z raasch
24! Bugfix: output of use_upstream_for_tke
25!
26! 82 2007-04-16 15:40:52Z raasch
27! Preprocessor strings for different linux clusters changed to "lc",
28! routine local_flush is used for buffer flushing
29!
30! 76 2007-03-29 00:58:32Z raasch
31! Output of netcdf_64bit_3d, particles-package is now part of the default code,
32! output of the loop optimization method, moisture renamed humidity,
33! output of subversion revision number
34!
35! 19 2007-02-23 04:53:48Z raasch
36! Output of scalar flux applied at top boundary
37!
38! RCS Log replace by Id keyword, revision history cleaned up
39!
40! Revision 1.63  2006/08/22 13:53:13  raasch
41! Output of dz_max
42!
43! Revision 1.1  1997/08/11 06:17:20  raasch
44! Initial revision
45!
46!
47! Description:
48! ------------
49! Writing a header with all important informations about the actual run.
50! This subroutine is called three times, two times at the beginning
51! (writing information on files RUN_CONTROL and HEADER) and one time at the
52! end of the run, then writing additional information about CPU-usage on file
53! header.
54!------------------------------------------------------------------------------!
55
56    USE arrays_3d
57    USE control_parameters
58    USE cloud_parameters
59    USE cpulog
60    USE dvrp_variables
61    USE grid_variables
62    USE indices
63    USE model_1d
64    USE particle_attributes
65    USE pegrid
66    USE spectrum
67
68    IMPLICIT NONE
69
70    CHARACTER (LEN=1)  ::  prec
71    CHARACTER (LEN=2)  ::  do2d_mode
72    CHARACTER (LEN=5)  ::  section_chr
73    CHARACTER (LEN=9)  ::  time_to_string
74    CHARACTER (LEN=10) ::  coor_chr, host_chr
75    CHARACTER (LEN=16) ::  begin_chr
76    CHARACTER (LEN=21) ::  ver_rev
77    CHARACTER (LEN=40) ::  output_format
78    CHARACTER (LEN=70) ::  char1, char2, coordinates, gradients, dopr_chr, &
79                           do2d_xy, do2d_xz, do2d_yz, do3d_chr, &
80                           run_classification, slices, temperatures, learde, &
81                           ugcomponent, vgcomponent
82    CHARACTER (LEN=85) ::  roben, runten
83
84    INTEGER ::  av, bh, blx, bly, bxl, bxr, byn, bys, i, ihost, io, j, l, ll
85    REAL    ::  cpuseconds_per_simulated_second
86
87!
88!-- Open the output file. At the end of the simulation, output is directed
89!-- to unit 19.
90    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
91         .NOT. simulated_time_at_begin /= simulated_time )  THEN
92       io = 15   !  header output on file RUN_CONTROL
93    ELSE
94       io = 19   !  header output on file HEADER
95    ENDIF
96    CALL check_open( io )
97
98!
99!-- At the end of the run, output file (HEADER) will be rewritten with
100!-- new informations
101    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
102
103!
104!-- Determine kind of model run
105    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
106       run_classification = '3D - restart run'
107    ELSE
108       IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
109          run_classification = '3D - run without 1D - prerun'
110       ELSEIF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 ) THEN
111          run_classification = '3D - run with 1D - prerun'
112       ELSE
113          PRINT*,'+++ header:  unknown action(s): ',initializing_actions
114       ENDIF
115    ENDIF
116    IF ( ocean )  THEN
117       run_classification = 'ocean - ' // run_classification
118    ELSE
119       run_classification = 'atmosphere - ' // run_classification
120    ENDIF
121
122!
123!-- Run-identification, date, time, host
124    host_chr = host(1:10)
125    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
126    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
127    IF ( coupling_mode /= 'uncoupled' )  WRITE ( io, 101 )  coupling_mode
128    WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr, &
129                       ADJUSTR( host_chr )
130#if defined( __parallel )
131    IF ( npex == -1  .AND.  pdims(2) /= 1 )  THEN
132       char1 = 'calculated'
133    ELSEIF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.  &
134               host(1:2) == 'lc' )  .AND.                          &
135             npex == -1  .AND.  pdims(2) == 1 )  THEN
136       char1 = 'forced'
137    ELSE
138       char1 = 'predefined'
139    ENDIF
140    IF ( threads_per_task == 1 )  THEN
141       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
142    ELSE
143       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
144                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
145    ENDIF
146    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
147           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
148         npex == -1  .AND.  pdims(2) == 1 )                      &
149    THEN
150       WRITE ( io, 106 )
151    ELSEIF ( pdims(2) == 1 )  THEN
152       WRITE ( io, 107 )  'x'
153    ELSEIF ( pdims(1) == 1 )  THEN
154       WRITE ( io, 107 )  'y'
155    ENDIF
156    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
157#endif
158    WRITE ( io, 99 )
159
160!
161!-- Numerical schemes
162    WRITE ( io, 110 )
163    IF ( psolver(1:7) == 'poisfft' )  THEN
164       WRITE ( io, 111 )  TRIM( fft_method )
165       IF ( psolver == 'poisfft_hybrid' )  WRITE ( io, 138 )
166    ELSEIF ( psolver == 'sor' )  THEN
167       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
168    ELSEIF ( psolver == 'multigrid' )  THEN
169       WRITE ( io, 135 )  cycle_mg, maximum_grid_level, ngsrb
170       IF ( mg_cycles == -1 )  THEN
171          WRITE ( io, 140 )  residual_limit
172       ELSE
173          WRITE ( io, 141 )  mg_cycles
174       ENDIF
175       IF ( mg_switch_to_pe0_level == 0 )  THEN
176          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
177                             nzt_mg(1)
178       ELSE
179          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
180                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
181                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
182                             nzt_mg(mg_switch_to_pe0_level),    &
183                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
184                             nzt_mg(1)
185       ENDIF
186    ENDIF
187    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
188    THEN
189       WRITE ( io, 142 )
190    ENDIF
191
192    IF ( momentum_advec == 'pw-scheme' )  THEN
193       WRITE ( io, 113 )
194    ELSE
195       WRITE ( io, 114 )
196       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
197       IF ( overshoot_limit_u /= 0.0  .OR.  overshoot_limit_v /= 0.0  .OR. &
198            overshoot_limit_w /= 0.0 )  THEN
199          WRITE ( io, 127 )  overshoot_limit_u, overshoot_limit_v, &
200                             overshoot_limit_w
201       ENDIF
202       IF ( ups_limit_u /= 0.0  .OR.  ups_limit_v /= 0.0  .OR. &
203            ups_limit_w /= 0.0 )                               &
204       THEN
205          WRITE ( io, 125 )  ups_limit_u, ups_limit_v, ups_limit_w
206       ENDIF
207       IF ( long_filter_factor /= 0.0 )  WRITE ( io, 115 )  long_filter_factor
208    ENDIF
209    IF ( scalar_advec == 'pw-scheme' )  THEN
210       WRITE ( io, 116 )
211    ELSEIF ( scalar_advec == 'ups-scheme' )  THEN
212       WRITE ( io, 117 )
213       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
214       IF ( overshoot_limit_e /= 0.0  .OR.  overshoot_limit_pt /= 0.0 )  THEN
215          WRITE ( io, 128 )  overshoot_limit_e, overshoot_limit_pt
216       ENDIF
217       IF ( ups_limit_e /= 0.0  .OR.  ups_limit_pt /= 0.0 )  THEN
218          WRITE ( io, 126 )  ups_limit_e, ups_limit_pt
219       ENDIF
220    ELSE
221       WRITE ( io, 118 )
222    ENDIF
223
224    WRITE ( io, 139 )  TRIM( loop_optimization )
225
226    IF ( galilei_transformation )  THEN
227       IF ( use_ug_for_galilei_tr )  THEN
228          char1 = 'geostrophic wind'
229       ELSE
230          char1 = 'mean wind in model domain'
231       ENDIF
232       IF ( simulated_time_at_begin == simulated_time )  THEN
233          char2 = 'at the start of the run'
234       ELSE
235          char2 = 'at the end of the run'
236       ENDIF
237       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ), &
238                          advected_distance_x/1000.0, advected_distance_y/1000.0
239    ENDIF
240    IF ( timestep_scheme == 'leapfrog' )  THEN
241       WRITE ( io, 120 )
242    ELSEIF ( timestep_scheme == 'leapfrog+euler' )  THEN
243       WRITE ( io, 121 )
244    ELSE
245       WRITE ( io, 122 )  timestep_scheme
246    ENDIF
247    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
248    IF ( rayleigh_damping_factor /= 0.0 )  THEN
249       IF ( .NOT. ocean )  THEN
250          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
251               rayleigh_damping_factor
252       ELSE
253          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
254               rayleigh_damping_factor
255       ENDIF
256    ENDIF
257    IF ( humidity )  THEN
258       IF ( .NOT. cloud_physics )  THEN
259          WRITE ( io, 129 )
260       ELSE
261          WRITE ( io, 130 )
262          WRITE ( io, 131 )
263          IF ( radiation )      WRITE ( io, 132 )
264          IF ( precipitation )  WRITE ( io, 133 )
265       ENDIF
266    ENDIF
267    IF ( passive_scalar )  WRITE ( io, 134 )
268    IF ( conserve_volume_flow )  WRITE ( io, 150 )
269    WRITE ( io, 99 )
270
271!
272!-- Runtime and timestep informations
273    WRITE ( io, 200 )
274    IF ( .NOT. dt_fixed )  THEN
275       WRITE ( io, 201 )  dt_max, cfl_factor
276    ELSE
277       WRITE ( io, 202 )  dt
278    ENDIF
279    WRITE ( io, 203 )  simulated_time_at_begin, end_time
280
281    IF ( time_restart /= 9999999.9  .AND. &
282         simulated_time_at_begin == simulated_time )  THEN
283       IF ( dt_restart == 9999999.9 )  THEN
284          WRITE ( io, 204 )  ' Restart at:       ',time_restart
285       ELSE
286          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
287       ENDIF
288    ENDIF
289
290    IF ( simulated_time_at_begin /= simulated_time )  THEN
291       i = MAX ( log_point_s(10)%counts, 1 )
292       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0 )  THEN
293          cpuseconds_per_simulated_second = 0.0
294       ELSE
295          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
296                                            ( simulated_time -    &
297                                              simulated_time_at_begin )
298       ENDIF
299       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum, &
300                          log_point_s(10)%sum / REAL( i ),     &
301                          cpuseconds_per_simulated_second
302       IF ( time_restart /= 9999999.9  .AND.  time_restart < end_time )  THEN
303          IF ( dt_restart == 9999999.9 )  THEN
304             WRITE ( io, 204 )  ' Next restart at:  ',time_restart
305          ELSE
306             WRITE ( io, 205 )  ' Next restart at:  ',time_restart, dt_restart
307          ENDIF
308       ENDIF
309    ENDIF
310
311!
312!-- Computational grid
313    IF ( .NOT. ocean )  THEN
314       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
315       IF ( dz_stretch_level_index < nzt+1 )  THEN
316          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
317                             dz_stretch_factor, dz_max
318       ENDIF
319    ELSE
320       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
321       IF ( dz_stretch_level_index > 0 )  THEN
322          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
323                             dz_stretch_factor, dz_max
324       ENDIF
325    ENDIF
326    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
327                       MIN( nnz+2, nzt+2 )
328    IF ( numprocs > 1 )  THEN
329       IF ( nxa == nx  .AND.  nya == ny  .AND.  nza == nz )  THEN
330          WRITE ( io, 255 )
331       ELSE
332          WRITE ( io, 256 )  nnx-(nxa-nx), nny-(nya-ny), nzt+2
333       ENDIF
334    ENDIF
335    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
336
337!
338!-- Topography
339    WRITE ( io, 270 )  topography
340    SELECT CASE ( TRIM( topography ) )
341
342       CASE ( 'flat' )
343          ! no actions necessary
344
345       CASE ( 'single_building' )
346          blx = INT( building_length_x / dx )
347          bly = INT( building_length_y / dy )
348          bh  = INT( building_height / dz )
349
350          IF ( building_wall_left == 9999999.9 )  THEN
351             building_wall_left = ( nx + 1 - blx ) / 2 * dx
352          ENDIF
353          bxl = INT ( building_wall_left / dx + 0.5 )
354          bxr = bxl + blx
355
356          IF ( building_wall_south == 9999999.9 )  THEN
357             building_wall_south = ( ny + 1 - bly ) / 2 * dy
358          ENDIF
359          bys = INT ( building_wall_south / dy + 0.5 )
360          byn = bys + bly
361
362          WRITE ( io, 271 )  building_length_x, building_length_y, &
363                             building_height, bxl, bxr, bys, byn
364
365    END SELECT
366
367    IF ( plant_canopy ) THEN
368
369       WRITE ( io, 280 ) canopy_mode, pch_index, drag_coefficient
370
371!
372!--    Leaf area density profile
373!--    Building output strings, starting with surface value
374       WRITE ( learde, '(F6.2)' )  lad_surface
375       gradients = '------'
376       slices = '     0'
377       coordinates = '   0.0'
378       i = 1
379       DO  WHILE ( lad_vertical_gradient_level_ind(i) /= -9999 )
380
381          WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
382          learde = TRIM( learde ) // ' ' // TRIM( coor_chr )
383
384          WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
385          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
386
387          WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
388          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
389
390          WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
391          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
392
393          i = i + 1
394       ENDDO
395
396       WRITE ( io, 281 )  TRIM( coordinates ), TRIM( learde ), &
397                          TRIM( gradients ), TRIM( slices )
398
399    ENDIF
400
401!
402!-- Boundary conditions
403    IF ( ibc_p_b == 0 )  THEN
404       runten = 'p(0)     = 0      |'
405    ELSEIF ( ibc_p_b == 1 )  THEN
406       runten = 'p(0)     = p(1)   |'
407    ELSE
408       runten = 'p(0)     = p(1) +R|'
409    ENDIF
410    IF ( ibc_p_t == 0 )  THEN
411       roben  = 'p(nzt+1) = 0      |'
412    ELSE
413       roben  = 'p(nzt+1) = p(nzt) |'
414    ENDIF
415
416    IF ( ibc_uv_b == 0 )  THEN
417       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
418    ELSE
419       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
420    ENDIF
421    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
422       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
423    ELSEIF ( ibc_uv_t == 0 )  THEN
424       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
425    ELSE
426       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
427    ENDIF
428
429    IF ( ibc_pt_b == 0 )  THEN
430       runten = TRIM( runten ) // ' pt(0)   = pt_surface'
431    ELSEIF ( ibc_pt_b == 1 )  THEN
432       runten = TRIM( runten ) // ' pt(0)   = pt(1)'
433    ELSEIF ( ibc_pt_b == 2 )  THEN
434       runten = TRIM( runten ) // ' pt(0) = from coupled model'
435    ENDIF
436    IF ( ibc_pt_t == 0 )  THEN
437       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
438    ELSEIF( ibc_pt_t == 1 )  THEN
439       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
440    ELSEIF( ibc_pt_t == 2 )  THEN
441       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
442    ENDIF
443
444    WRITE ( io, 300 )  runten, roben
445
446    IF ( .NOT. constant_diffusion )  THEN
447       IF ( ibc_e_b == 1 )  THEN
448          runten = 'e(0)     = e(1)'
449       ELSE
450          runten = 'e(0)     = e(1) = (u*/0.1)**2'
451       ENDIF
452       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
453
454       WRITE ( io, 301 )  'e', runten, roben       
455
456    ENDIF
457
458    IF ( ocean )  THEN
459       runten = 'sa(0)    = sa(1)'
460       IF ( ibc_sa_t == 0 )  THEN
461          roben =  'sa(nzt+1) = sa_surface'
462       ELSE
463          roben =  'sa(nzt+1) = sa(nzt)'
464       ENDIF
465       WRITE ( io, 301 ) 'sa', runten, roben
466    ENDIF
467
468    IF ( humidity )  THEN
469       IF ( ibc_q_b == 0 )  THEN
470          runten = 'q(0)     = q_surface'
471       ELSE
472          runten = 'q(0)     = q(1)'
473       ENDIF
474       IF ( ibc_q_t == 0 )  THEN
475          roben =  'q(nzt)   = q_top'
476       ELSE
477          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
478       ENDIF
479       WRITE ( io, 301 ) 'q', runten, roben
480    ENDIF
481
482    IF ( passive_scalar )  THEN
483       IF ( ibc_q_b == 0 )  THEN
484          runten = 's(0)     = s_surface'
485       ELSE
486          runten = 's(0)     = s(1)'
487       ENDIF
488       IF ( ibc_q_t == 0 )  THEN
489          roben =  's(nzt)   = s_top'
490       ELSE
491          roben =  's(nzt)   = s(nzt-1) + ds/dz'
492       ENDIF
493       WRITE ( io, 301 ) 's', runten, roben
494    ENDIF
495
496    IF ( use_surface_fluxes )  THEN
497       WRITE ( io, 303 )
498       IF ( constant_heatflux )  THEN
499          WRITE ( io, 306 )  surface_heatflux
500          IF ( random_heatflux )  WRITE ( io, 307 )
501       ENDIF
502       IF ( humidity  .AND.  constant_waterflux )  THEN
503          WRITE ( io, 311 ) surface_waterflux
504       ENDIF
505       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
506          WRITE ( io, 313 ) surface_waterflux
507       ENDIF
508    ENDIF
509
510    IF ( use_top_fluxes )  THEN
511       WRITE ( io, 304 )
512       IF ( coupling_mode == 'uncoupled' )  THEN
513          WRITE ( io, 319 )  top_momentumflux_u, top_momentumflux_v
514          IF ( constant_top_heatflux )  THEN
515             WRITE ( io, 306 )  top_heatflux
516          ENDIF
517       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
518          WRITE ( io, 316 )
519       ENDIF
520       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
521          WRITE ( io, 309 )  top_salinityflux
522       ENDIF
523       IF ( humidity  .OR.  passive_scalar )  THEN
524          WRITE ( io, 315 )
525       ENDIF
526    ENDIF
527
528    IF ( prandtl_layer )  THEN
529       WRITE ( io, 305 )  0.5 * (zu(1)-zu(0)), roughness_length, kappa, &
530                          rif_min, rif_max
531       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
532       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
533          WRITE ( io, 312 )
534       ENDIF
535       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
536          WRITE ( io, 314 )
537       ENDIF
538    ELSE
539       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
540          WRITE ( io, 310 )  rif_min, rif_max
541       ENDIF
542    ENDIF
543
544    WRITE ( io, 317 )  bc_lr, bc_ns
545    IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
546       WRITE ( io, 318 )  outflow_damping_width, km_damp_max
547    ENDIF
548
549!
550!-- Listing of 1D-profiles
551    WRITE ( io, 320 )  dt_dopr_listing
552    IF ( averaging_interval_pr /= 0.0 )  THEN
553       WRITE ( io, 321 )  averaging_interval_pr, dt_averaging_input_pr
554    ENDIF
555
556!
557!-- DATA output
558    WRITE ( io, 330 )
559    IF ( averaging_interval_pr /= 0.0 )  THEN
560       WRITE ( io, 321 )  averaging_interval_pr, dt_averaging_input_pr
561    ENDIF
562
563!
564!-- 1D-profiles
565    dopr_chr = 'Profile:'
566    IF ( dopr_n /= 0 )  THEN
567       WRITE ( io, 331 )
568
569       output_format = ''
570       IF ( netcdf_output )  THEN
571          IF ( netcdf_64bit )  THEN
572             output_format = 'netcdf (64 bit offset)'
573          ELSE
574             output_format = 'netcdf'
575          ENDIF
576       ENDIF
577       IF ( profil_output )  THEN
578          IF ( netcdf_output )  THEN
579             output_format = TRIM( output_format ) // ' and profil'
580          ELSE
581             output_format = 'profil'
582          ENDIF
583       ENDIF
584       WRITE ( io, 345 )  output_format
585
586       DO  i = 1, dopr_n
587          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
588          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
589             WRITE ( io, 332 )  dopr_chr
590             dopr_chr = '       :'
591          ENDIF
592       ENDDO
593
594       IF ( dopr_chr /= '' )  THEN
595          WRITE ( io, 332 )  dopr_chr
596       ENDIF
597       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
598       IF ( skip_time_dopr /= 0.0 )  WRITE ( io, 339 )  skip_time_dopr
599    ENDIF
600
601!
602!-- 2D-arrays
603    DO  av = 0, 1
604
605       i = 1
606       do2d_xy = ''
607       do2d_xz = ''
608       do2d_yz = ''
609       DO  WHILE ( do2d(av,i) /= ' ' )
610
611          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
612          do2d_mode = do2d(av,i)(l-1:l)
613
614          SELECT CASE ( do2d_mode )
615             CASE ( 'xy' )
616                ll = LEN_TRIM( do2d_xy )
617                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
618             CASE ( 'xz' )
619                ll = LEN_TRIM( do2d_xz )
620                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
621             CASE ( 'yz' )
622                ll = LEN_TRIM( do2d_yz )
623                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
624          END SELECT
625
626          i = i + 1
627
628       ENDDO
629
630       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
631              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
632              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) )  .AND. &
633            ( netcdf_output  .OR.  iso2d_output ) )  THEN
634
635          IF (  av == 0 )  THEN
636             WRITE ( io, 334 )  ''
637          ELSE
638             WRITE ( io, 334 )  '(time-averaged)'
639          ENDIF
640
641          IF ( do2d_at_begin )  THEN
642             begin_chr = 'and at the start'
643          ELSE
644             begin_chr = ''
645          ENDIF
646
647          output_format = ''
648          IF ( netcdf_output )  THEN
649             IF ( netcdf_64bit )  THEN
650                output_format = 'netcdf (64 bit offset)'
651             ELSE
652                output_format = 'netcdf'
653             ENDIF
654          ENDIF
655          IF ( iso2d_output )  THEN
656             IF ( netcdf_output )  THEN
657                output_format = TRIM( output_format ) // ' and iso2d'
658             ELSE
659                output_format = 'iso2d'
660             ENDIF
661          ENDIF
662          WRITE ( io, 345 )  output_format
663
664          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
665             i = 1
666             slices = '/'
667             coordinates = '/'
668!
669!--          Building strings with index and coordinate informations of the
670!--          slices
671             DO  WHILE ( section(i,1) /= -9999 )
672
673                WRITE (section_chr,'(I5)')  section(i,1)
674                section_chr = ADJUSTL( section_chr )
675                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
676
677                WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
678                coor_chr = ADJUSTL( coor_chr )
679                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
680
681                i = i + 1
682             ENDDO
683             IF ( av == 0 )  THEN
684                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
685                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
686                                   TRIM( coordinates )
687                IF ( skip_time_do2d_xy /= 0.0 )  THEN
688                   WRITE ( io, 339 )  skip_time_do2d_xy
689                ENDIF
690             ELSE
691                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
692                                   TRIM( begin_chr ), averaging_interval, &
693                                   dt_averaging_input, 'k', TRIM( slices ), &
694                                   TRIM( coordinates )
695                IF ( skip_time_data_output_av /= 0.0 )  THEN
696                   WRITE ( io, 339 )  skip_time_data_output_av
697                ENDIF
698             ENDIF
699
700          ENDIF
701
702          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
703             i = 1
704             slices = '/'
705             coordinates = '/'
706!
707!--          Building strings with index and coordinate informations of the
708!--          slices
709             DO  WHILE ( section(i,2) /= -9999 )
710
711                WRITE (section_chr,'(I5)')  section(i,2)
712                section_chr = ADJUSTL( section_chr )
713                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
714
715                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
716                coor_chr = ADJUSTL( coor_chr )
717                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
718
719                i = i + 1
720             ENDDO
721             IF ( av == 0 )  THEN
722                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
723                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
724                                   TRIM( coordinates )
725                IF ( skip_time_do2d_xz /= 0.0 )  THEN
726                   WRITE ( io, 339 )  skip_time_do2d_xz
727                ENDIF
728             ELSE
729                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
730                                   TRIM( begin_chr ), averaging_interval, &
731                                   dt_averaging_input, 'j', TRIM( slices ), &
732                                   TRIM( coordinates )
733                IF ( skip_time_data_output_av /= 0.0 )  THEN
734                   WRITE ( io, 339 )  skip_time_data_output_av
735                ENDIF
736             ENDIF
737          ENDIF
738
739          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
740             i = 1
741             slices = '/'
742             coordinates = '/'
743!
744!--          Building strings with index and coordinate informations of the
745!--          slices
746             DO  WHILE ( section(i,3) /= -9999 )
747
748                WRITE (section_chr,'(I5)')  section(i,3)
749                section_chr = ADJUSTL( section_chr )
750                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
751
752                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
753                coor_chr = ADJUSTL( coor_chr )
754                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
755
756                i = i + 1
757             ENDDO
758             IF ( av == 0 )  THEN
759                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
760                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
761                                   TRIM( coordinates )
762                IF ( skip_time_do2d_yz /= 0.0 )  THEN
763                   WRITE ( io, 339 )  skip_time_do2d_yz
764                ENDIF
765             ELSE
766                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
767                                   TRIM( begin_chr ), averaging_interval, &
768                                   dt_averaging_input, 'i', TRIM( slices ), &
769                                   TRIM( coordinates )
770                IF ( skip_time_data_output_av /= 0.0 )  THEN
771                   WRITE ( io, 339 )  skip_time_data_output_av
772                ENDIF
773             ENDIF
774          ENDIF
775
776       ENDIF
777
778    ENDDO
779
780!
781!-- 3d-arrays
782    DO  av = 0, 1
783
784       i = 1
785       do3d_chr = ''
786       DO  WHILE ( do3d(av,i) /= ' ' )
787
788          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
789          i = i + 1
790
791       ENDDO
792
793       IF ( do3d_chr /= '' )  THEN
794          IF ( av == 0 )  THEN
795             WRITE ( io, 336 )  ''
796          ELSE
797             WRITE ( io, 336 )  '(time-averaged)'
798          ENDIF
799
800          output_format = ''
801          IF ( netcdf_output )  THEN
802             IF ( netcdf_64bit .AND. netcdf_64bit_3d )  THEN
803                output_format = 'netcdf (64 bit offset)'
804             ELSE
805                output_format = 'netcdf'
806             ENDIF
807          ENDIF
808          IF ( avs_output )  THEN
809             IF ( netcdf_output )  THEN
810                output_format = TRIM( output_format ) // ' and avs'
811             ELSE
812                output_format = 'avs'
813             ENDIF
814          ENDIF
815          WRITE ( io, 345 )  output_format
816
817          IF ( do3d_at_begin )  THEN
818             begin_chr = 'and at the start'
819          ELSE
820             begin_chr = ''
821          ENDIF
822          IF ( av == 0 )  THEN
823             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
824                                zu(nz_do3d), nz_do3d
825          ELSE
826             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
827                                TRIM( begin_chr ), averaging_interval, &
828                                dt_averaging_input, zu(nz_do3d), nz_do3d
829          ENDIF
830
831          IF ( do3d_compress )  THEN
832             do3d_chr = ''
833             i = 1
834             DO WHILE ( do3d(av,i) /= ' ' )
835
836                SELECT CASE ( do3d(av,i) )
837                   CASE ( 'u' )
838                      j = 1
839                   CASE ( 'v' )
840                      j = 2
841                   CASE ( 'w' )
842                      j = 3
843                   CASE ( 'p' )
844                      j = 4
845                   CASE ( 'pt' )
846                      j = 5
847                END SELECT
848                WRITE ( prec, '(I1)' )  plot_3d_precision(j)%precision
849                do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // &
850                           ':' // prec // ','
851                i = i + 1
852
853             ENDDO
854             WRITE ( io, 338 )  do3d_chr
855
856          ENDIF
857
858          IF ( av == 0 )  THEN
859             IF ( skip_time_do3d /= 0.0 )  THEN
860                WRITE ( io, 339 )  skip_time_do3d
861             ENDIF
862          ELSE
863             IF ( skip_time_data_output_av /= 0.0 )  THEN
864                WRITE ( io, 339 )  skip_time_data_output_av
865             ENDIF
866          ENDIF
867
868       ENDIF
869
870    ENDDO
871
872!
873!-- Timeseries
874    IF ( dt_dots /= 9999999.9 )  THEN
875       WRITE ( io, 340 )
876
877       output_format = ''
878       IF ( netcdf_output )  THEN
879          IF ( netcdf_64bit )  THEN
880             output_format = 'netcdf (64 bit offset)'
881          ELSE
882             output_format = 'netcdf'
883          ENDIF
884       ENDIF
885       IF ( profil_output )  THEN
886          IF ( netcdf_output )  THEN
887             output_format = TRIM( output_format ) // ' and profil'
888          ELSE
889             output_format = 'profil'
890          ENDIF
891       ENDIF
892       WRITE ( io, 345 )  output_format
893       WRITE ( io, 341 )  dt_dots
894    ENDIF
895
896#if defined( __dvrp_graphics )
897!
898!-- Dvrp-output
899    IF ( dt_dvrp /= 9999999.9 )  THEN
900       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
901                          TRIM( dvrp_username ), TRIM( dvrp_directory )
902       i = 1
903       l = 0
904       DO WHILE ( mode_dvrp(i) /= ' ' )
905          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
906             READ ( mode_dvrp(i), '(10X,I2)' )  j
907             l = l + 1
908             IF ( do3d(0,j) /= ' ' )  THEN
909                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l)
910             ENDIF
911          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
912             READ ( mode_dvrp(i), '(6X,I2)' )  j
913             IF ( do2d(0,j) /= ' ' )  WRITE ( io, 362 )  TRIM( do2d(0,j) )
914          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
915             WRITE ( io, 363 )
916          ENDIF
917          i = i + 1
918       ENDDO
919    ENDIF
920#endif
921
922#if defined( __spectra )
923!
924!-- Spectra output
925    IF ( dt_dosp /= 9999999.9 ) THEN
926       WRITE ( io, 370 )
927
928       output_format = ''
929       IF ( netcdf_output )  THEN
930          IF ( netcdf_64bit )  THEN
931             output_format = 'netcdf (64 bit offset)'
932          ELSE
933             output_format = 'netcdf'
934          ENDIF
935       ENDIF
936       IF ( profil_output )  THEN
937          IF ( netcdf_output )  THEN
938             output_format = TRIM( output_format ) // ' and profil'
939          ELSE
940             output_format = 'profil'
941          ENDIF
942       ENDIF
943       WRITE ( io, 345 )  output_format
944       WRITE ( io, 371 )  dt_dosp
945       IF ( skip_time_dosp /= 0.0 )  WRITE ( io, 339 )  skip_time_dosp
946       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
947                          ( spectra_direction(i), i = 1,10 ),  &
948                          ( comp_spectra_level(i), i = 1,10 ), &
949                          ( plot_spectra_level(i), i = 1,10 ), &
950                          averaging_interval_sp, dt_averaging_input_pr
951    ENDIF
952#endif
953
954    WRITE ( io, 99 )
955
956!
957!-- Physical quantities
958    WRITE ( io, 400 )
959
960!
961!-- Geostrophic parameters
962    WRITE ( io, 410 )  omega, phi, f, fs
963
964!
965!-- Other quantities
966    WRITE ( io, 411 )  g
967    IF ( use_reference )  THEN
968       IF ( ocean )  THEN
969          WRITE ( io, 412 )  prho_reference
970       ELSE
971          WRITE ( io, 413 )  pt_reference
972       ENDIF
973    ENDIF
974
975!
976!-- Cloud physics parameters
977    IF ( cloud_physics ) THEN
978       WRITE ( io, 415 )
979       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
980    ENDIF
981
982!-- Profile of the geostrophic wind (component ug)
983!-- Building output strings
984    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
985    gradients = '------'
986    slices = '     0'
987    coordinates = '   0.0'
988    i = 1
989    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
990     
991       WRITE (coor_chr,'(F6.2,4X)')  ug(ug_vertical_gradient_level_ind(i))
992       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
993
994       WRITE (coor_chr,'(F6.2,4X)')  ug_vertical_gradient(i)
995       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
996
997       WRITE (coor_chr,'(I6,4X)')  ug_vertical_gradient_level_ind(i)
998       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
999
1000       WRITE (coor_chr,'(F6.1,4X)')  ug_vertical_gradient_level(i)
1001       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1002
1003       i = i + 1
1004    ENDDO
1005
1006    WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
1007                       TRIM( gradients ), TRIM( slices )
1008
1009!-- Profile of the geostrophic wind (component vg)
1010!-- Building output strings
1011    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
1012    gradients = '------'
1013    slices = '     0'
1014    coordinates = '   0.0'
1015    i = 1
1016    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
1017
1018       WRITE (coor_chr,'(F6.2,4X)')  vg(vg_vertical_gradient_level_ind(i))
1019       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
1020
1021       WRITE (coor_chr,'(F6.2,4X)')  vg_vertical_gradient(i)
1022       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1023
1024       WRITE (coor_chr,'(I6,4X)')  vg_vertical_gradient_level_ind(i)
1025       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1026
1027       WRITE (coor_chr,'(F6.1,4X)')  vg_vertical_gradient_level(i)
1028       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1029
1030       i = i + 1 
1031    ENDDO
1032
1033    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
1034                       TRIM( gradients ), TRIM( slices )
1035
1036!
1037!-- Initial temperature profile
1038!-- Building output strings, starting with surface temperature
1039    WRITE ( temperatures, '(F6.2)' )  pt_surface
1040    gradients = '------'
1041    slices = '     0'
1042    coordinates = '   0.0'
1043    i = 1
1044    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1045
1046       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1047       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1048
1049       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1050       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1051
1052       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1053       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1054
1055       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1056       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1057
1058       i = i + 1
1059    ENDDO
1060
1061    WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1062                       TRIM( gradients ), TRIM( slices )
1063
1064!
1065!-- Initial humidity profile
1066!-- Building output strings, starting with surface humidity
1067    IF ( humidity  .OR.  passive_scalar )  THEN
1068       WRITE ( temperatures, '(E8.1)' )  q_surface
1069       gradients = '--------'
1070       slices = '       0'
1071       coordinates = '     0.0'
1072       i = 1
1073       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1074         
1075          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1076          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1077
1078          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1079          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1080         
1081          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1082          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1083         
1084          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1085          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1086
1087          i = i + 1
1088       ENDDO
1089
1090       IF ( humidity )  THEN
1091          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1092                             TRIM( gradients ), TRIM( slices )
1093       ELSE
1094          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1095                             TRIM( gradients ), TRIM( slices )
1096       ENDIF
1097    ENDIF
1098
1099!
1100!-- Initial salinity profile
1101!-- Building output strings, starting with surface salinity
1102    IF ( ocean )  THEN
1103       WRITE ( temperatures, '(F6.2)' )  sa_surface
1104       gradients = '------'
1105       slices = '     0'
1106       coordinates = '   0.0'
1107       i = 1
1108       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1109
1110          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1111          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1112
1113          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1114          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1115
1116          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1117          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1118
1119          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1120          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1121
1122          i = i + 1
1123       ENDDO
1124
1125       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1126                          TRIM( gradients ), TRIM( slices )
1127    ENDIF
1128
1129!
1130!-- LES / turbulence parameters
1131    WRITE ( io, 450 )
1132
1133!--
1134! ... LES-constants used must still be added here
1135!--
1136    IF ( constant_diffusion )  THEN
1137       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1138                          prandtl_number
1139    ENDIF
1140    IF ( .NOT. constant_diffusion)  THEN
1141       IF ( e_init > 0.0 )  WRITE ( io, 455 )  e_init
1142       IF ( e_min > 0.0 )  WRITE ( io, 454 )  e_min
1143       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1144       IF ( adjust_mixing_length  .AND.  prandtl_layer )  WRITE ( io, 452 )
1145    ENDIF
1146
1147!
1148!-- Special actions during the run
1149    WRITE ( io, 470 )
1150    IF ( create_disturbances )  THEN
1151       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1152                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1153                          zu(disturbance_level_ind_t), disturbance_level_ind_t
1154       IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
1155          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1156       ELSE
1157          WRITE ( io, 473 )  disturbance_energy_limit
1158       ENDIF
1159       WRITE ( io, 474 )  TRIM( random_generator )
1160    ENDIF
1161    IF ( pt_surface_initial_change /= 0.0 )  THEN
1162       WRITE ( io, 475 )  pt_surface_initial_change
1163    ENDIF
1164    IF ( humidity  .AND.  q_surface_initial_change /= 0.0 )  THEN
1165       WRITE ( io, 476 )  q_surface_initial_change       
1166    ENDIF
1167    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0 )  THEN
1168       WRITE ( io, 477 )  q_surface_initial_change       
1169    ENDIF
1170
1171    IF ( particle_advection )  THEN
1172!
1173!--    Particle attributes
1174       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1175                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
1176                          end_time_prel, dt_sort_particles
1177       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1178       IF ( random_start_position )  WRITE ( io, 481 )
1179       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1180       WRITE ( io, 495 )  total_number_of_particles
1181       IF ( .NOT. vertical_particle_advection )  WRITE ( io, 482 )
1182       IF ( maximum_number_of_tailpoints /= 0 )  THEN
1183          WRITE ( io, 483 )  maximum_number_of_tailpoints
1184          IF ( minimum_tailpoint_distance /= 0 )  THEN
1185             WRITE ( io, 484 )  total_number_of_tails,      &
1186                                minimum_tailpoint_distance, &
1187                                maximum_tailpoint_age
1188          ENDIF
1189       ENDIF
1190       IF ( dt_write_particle_data /= 9999999.9 )  THEN
1191          WRITE ( io, 485 )  dt_write_particle_data
1192          output_format = ''
1193          IF ( netcdf_output )  THEN
1194             IF ( netcdf_64bit )  THEN
1195                output_format = 'netcdf (64 bit offset) and binary'
1196             ELSE
1197                output_format = 'netcdf and binary'
1198             ENDIF
1199          ELSE
1200             output_format = 'binary'
1201          ENDIF
1202          WRITE ( io, 345 )  output_format
1203       ENDIF
1204       IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1205       IF ( write_particle_statistics )  WRITE ( io, 486 )
1206
1207       WRITE ( io, 487 )  number_of_particle_groups
1208
1209       DO  i = 1, number_of_particle_groups
1210          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9 )  THEN
1211             WRITE ( io, 490 )  i, 0.0
1212             WRITE ( io, 492 )
1213          ELSE
1214             WRITE ( io, 490 )  i, radius(i)
1215             IF ( density_ratio(i) /= 0.0 )  THEN
1216                WRITE ( io, 491 )  density_ratio(i)
1217             ELSE
1218                WRITE ( io, 492 )
1219             ENDIF
1220          ENDIF
1221          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1222                             pdx(i), pdy(i), pdz(i)
1223       ENDDO
1224
1225    ENDIF
1226
1227
1228!
1229!-- Parameters of 1D-model
1230    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1231       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1232                          mixing_length_1d, dissipation_1d
1233       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1234          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1235       ENDIF
1236    ENDIF
1237
1238!
1239!-- User-defined informations
1240    CALL user_header( io )
1241
1242    WRITE ( io, 99 )
1243
1244!
1245!-- Write buffer contents to disc immediately
1246    CALL local_flush( io )
1247
1248!
1249!-- Here the FORMATs start
1250
1251 99 FORMAT (1X,78('-'))
1252100 FORMAT (/1X,'*************************',11X,42('-')/        &
1253            1X,'* ',A,' *',11X,A/                               &
1254            1X,'*************************',11X,42('-'))
1255101 FORMAT (37X,'coupled run: ',A/ &
1256            37X,42('-'))
1257102 FORMAT (/' Date:            ',A8,11X,'Run:       ',A20/      &
1258            ' Time:            ',A8,11X,'Run-No.:   ',I2.2/     &
1259            ' Run on host:   ',A10)
1260#if defined( __parallel )
1261103 FORMAT (' Number of PEs:',7X,I4,11X,'Processor grid (x,y): (',I3,',',I3, &
1262              ')',1X,A)
1263104 FORMAT (' Number of PEs:',7X,I4,11X,'Tasks:',I4,'   threads per task:',I4/ &
1264              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
1265105 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1266106 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
1267            37X,'because the job is running on an SMP-cluster')
1268107 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
1269#endif
1270110 FORMAT (/' Numerical Schemes:'/ &
1271             ' -----------------'/)
1272111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1273112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1274            '     Iterations (initial/other): ',I3,'/',I3,'  omega = ',F5.3)
1275113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1276                  ' or Upstream')
1277114 FORMAT (' --> Momentum advection via Upstream-Spline-Scheme')
1278115 FORMAT ('     Tendencies are smoothed via Long-Filter with factor ',F5.3) 
1279116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1280                  ' or Upstream')
1281117 FORMAT (' --> Scalar advection via Upstream-Spline-Scheme')
1282118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1283119 FORMAT (' --> Galilei-Transform applied to horizontal advection', &
1284            '     Translation velocity = ',A/ &
1285            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1286120 FORMAT (' --> Time differencing scheme: leapfrog only (no euler in case', &
1287                  ' of timestep changes)')
1288121 FORMAT (' --> Time differencing scheme: leapfrog + euler in case of', &
1289                  ' timestep changes')
1290122 FORMAT (' --> Time differencing scheme: ',A)
1291123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
1292            '     maximum damping coefficient: ',F5.3, ' 1/s')
1293124 FORMAT ('     Spline-overshoots are being suppressed')
1294125 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1295                  ' of'/                                                       &
1296            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1297126 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1298                  ' of'/                                                       &
1299            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1300127 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1301            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1302128 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1303            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1304129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1305130 FORMAT (' --> Additional prognostic equation for the total water content')
1306131 FORMAT (' --> Parameterization of condensation processes via (0%-or100%)')
1307132 FORMAT (' --> Parameterization of long-wave radiation processes via'/ &
1308            '     effective emissivity scheme')
1309133 FORMAT (' --> Precipitation parameterization via Kessler-Scheme')
1310134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1311135 FORMAT (' --> Solve perturbation pressure via multigrid method (', &
1312                  A,'-cycle)'/ &
1313            '     number of grid levels:                   ',I2/ &
1314            '     Gauss-Seidel red/black iterations:       ',I2)
1315136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1316                  I3,')')
1317137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1318            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1319                  I3,')'/ &
1320            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1321                  I3,')')
1322138 FORMAT ('     Using hybrid version for 1d-domain-decomposition')
1323139 FORMAT (' --> Loop optimization method: ',A)
1324140 FORMAT ('     maximum residual allowed:                ',E10.3)
1325141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1326142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1327                  'step')
1328143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1329                  'kinetic energy')
1330150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
1331                  'conserved')
1332200 FORMAT (//' Run time and time step information:'/ &
1333             ' ----------------------------------'/)
1334201 FORMAT ( ' Timestep:          variable     maximum value: ',F6.3,' s', &
1335             '    CFL-factor: ',F4.2)
1336202 FORMAT ( ' Timestep:       dt = ',F6.3,' s'/)
1337203 FORMAT ( ' Start time:       ',F9.3,' s'/ &
1338             ' End time:         ',F9.3,' s')
1339204 FORMAT ( A,F9.3,' s')
1340205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
1341206 FORMAT (/' Time reached:     ',F9.3,' s'/ &
1342             ' CPU-time used:    ',F9.3,' s     per timestep:               ', &
1343               '  ',F9.3,' s'/                                                 &
1344             '                                   per second of simulated tim', &
1345               'e: ',F9.3,' s')
1346250 FORMAT (//' Computational grid and domain size:'/ &
1347              ' ----------------------------------'// &
1348              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1349              ' m    dz =    ',F7.3,' m'/ &
1350              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1351              ' m  z(u) = ',F10.3,' m'/)
1352252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
1353              ' factor: ',F5.3/ &
1354            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1355254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1356            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1357255 FORMAT (' Subdomains have equal size')
1358256 FORMAT (' Subdomains at the upper edges of the virtual processor grid ', &
1359              'have smaller sizes'/                                          &
1360            ' Size of smallest subdomain:    (  ',I4,',   ',I4,',   ',I4,')')
1361260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1362             ' degrees')
1363270 FORMAT (//' Topography informations:'/ &
1364              ' -----------------------'// &
1365              1X,'Topography: ',A)
1366271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1367              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1368                ' / ',I4)
1369280 FORMAT (//' Vegetation canopy (drag) model:'/ &
1370              ' ------------------------------'// &
1371              ' Canopy mode: ', A / &
1372              ' Canopy top: ',I4 / &
1373              ' Leaf drag coefficient: ',F6.2 /)
1374281 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
1375              ' Height:              ',A,'  m'/ &
1376              ' Leaf area density:   ',A,'  m**2/m**3'/ &
1377              ' Gradient:            ',A,'  m**2/m**4'/ &
1378              ' Gridpoint:           ',A)
1379               
1380300 FORMAT (//' Boundary conditions:'/ &
1381             ' -------------------'// &
1382             '                     p                    uv             ', &
1383             '                   pt'// &
1384             ' B. bound.: ',A/ &
1385             ' T. bound.: ',A)
1386301 FORMAT (/'                     ',A// &
1387             ' B. bound.: ',A/ &
1388             ' T. bound.: ',A)
1389303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1390304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1391305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
1392               'computational u,v-level:'// &
1393             '       zp = ',F6.2,' m   z0 = ',F6.4,' m   kappa = ',F4.2/ &
1394             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1395306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
1396307 FORMAT ('       Heatflux has a random normal distribution')
1397308 FORMAT ('       Predefined surface temperature')
1398309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
1399310 FORMAT (//'    1D-Model:'// &
1400             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1401311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1402312 FORMAT ('       Predefined surface humidity')
1403313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1404314 FORMAT ('       Predefined scalar value at the surface')
1405315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
1406316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1407                    'atmosphere model')
1408317 FORMAT (//' Lateral boundaries:'/ &
1409            '       left/right:  ',A/    &
1410            '       north/south: ',A)
1411318 FORMAT (/'       outflow damping layer width: ',I3,' gridpoints with km_', &
1412                    'max =',F5.1,' m**2/s')
1413319 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
1414            '                                          v: ',F9.6,' m**2/s**2')
1415320 FORMAT (//' List output:'/ &
1416             ' -----------'//  &
1417            '    1D-Profiles:'/    &
1418            '       Output every             ',F8.2,' s')
1419321 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
1420            '       Averaging input every    ',F8.2,' s')
1421330 FORMAT (//' Data output:'/ &
1422             ' -----------'/)
1423331 FORMAT (/'    1D-Profiles:')
1424332 FORMAT (/'       ',A)
1425333 FORMAT ('       Output every             ',F8.2,' s',/ &
1426            '       Time averaged over       ',F8.2,' s'/ &
1427            '       Averaging input every    ',F8.2,' s')
1428334 FORMAT (/'    2D-Arrays',A,':')
1429335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1430            '       Output every             ',F8.2,' s  ',A/ &
1431            '       Cross sections at ',A1,' = ',A/ &
1432            '       scalar-coordinates:   ',A,' m'/)
1433336 FORMAT (/'    3D-Arrays',A,':')
1434337 FORMAT (/'       Arrays: ',A/ &
1435            '       Output every             ',F8.2,' s  ',A/ &
1436            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1437338 FORMAT ('       Compressed data output'/ &
1438            '       Decimal precision: ',A/)
1439339 FORMAT ('       No output during initial ',F8.2,' s')
1440340 FORMAT (/'    Time series:')
1441341 FORMAT ('       Output every             ',F8.2,' s'/)
1442342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1443            '       Output every             ',F8.2,' s  ',A/ &
1444            '       Time averaged over       ',F8.2,' s'/ &
1445            '       Averaging input every    ',F8.2,' s'/ &
1446            '       Cross sections at ',A1,' = ',A/ &
1447            '       scalar-coordinates:   ',A,' m'/)
1448343 FORMAT (/'       Arrays: ',A/ &
1449            '       Output every             ',F8.2,' s  ',A/ &
1450            '       Time averaged over       ',F8.2,' s'/ &
1451            '       Averaging input every    ',F8.2,' s'/ &
1452            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1453345 FORMAT ('       Output format: ',A/)
1454#if defined( __dvrp_graphics )
1455360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
1456            '       Output every      ',F7.1,' s'/ &
1457            '       Output mode:      ',A/ &
1458            '       Host / User:      ',A,' / ',A/ &
1459            '       Directory:        ',A// &
1460            '       The sequence contains:')
1461361 FORMAT ('       Isosurface of ',A,'  Threshold value: ', E12.3)
1462362 FORMAT ('       Sectional plane ',A)
1463363 FORMAT ('       Particles')
1464#endif
1465#if defined( __spectra )
1466370 FORMAT ('    Spectra:')
1467371 FORMAT ('       Output every ',F7.1,' s'/)
1468372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
1469            '       Directions: ', 10(A5,',')/                         &
1470            '       height levels  k = ', 9(I3,','),I3,'.'/            &
1471            '       height levels selected for standard plot:'/        &
1472            '                      k = ', 9(I3,','),I3,'.'/            &
1473            '       Time averaged over ', F7.1, ' s,' /                &
1474            '       Profiles for the time averaging are taken every ', &
1475                    F6.1,' s')
1476#endif
1477400 FORMAT (//' Physical quantities:'/ &
1478              ' -------------------'/)
1479410 FORMAT ('    Angular velocity    :   omega = ',E9.3,' rad/s'/  &
1480            '    Geograph. latitude  :   phi   = ',F4.1,' degr'/   &
1481            '    Coriolis parameter  :   f     = ',F9.6,' 1/s'/    &
1482            '                            f*    = ',F9.6,' 1/s')
1483411 FORMAT (/'    Gravity             :   g     = ',F4.1,' m/s**2')
1484412 FORMAT (/'    Reference density in buoyancy terms: ',F8.3,' kg/m**3')
1485413 FORMAT (/'    Reference temperature in buoyancy terms: ',F8.4,' K')
1486415 FORMAT (/'    Cloud physics parameters:'/ &
1487             '    ------------------------'/)
1488416 FORMAT ('        Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
1489            '        Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
1490            '        Density of air     :   rho_0 = ',F5.3,' kg/m**3'/  &
1491            '        Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
1492            '        Vapourization heat :   L_v   = ',E8.2,' J/kg')
1493420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1494            '       Height:        ',A,'  m'/ &
1495            '       Temperature:   ',A,'  K'/ &
1496            '       Gradient:      ',A,'  K/100m'/ &
1497            '       Gridpoint:     ',A)
1498421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1499            '       Height:      ',A,'  m'/ &
1500            '       Humidity:    ',A,'  kg/kg'/ &
1501            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1502            '       Gridpoint:   ',A)
1503422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1504            '       Height:                  ',A,'  m'/ &
1505            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1506            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1507            '       Gridpoint:               ',A)
1508423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1509            '       Height:      ',A,'  m'/ &
1510            '       ug:          ',A,'  m/s'/ &
1511            '       Gradient:    ',A,'  1/100s'/ &
1512            '       Gridpoint:   ',A)
1513424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1514            '       Height:      ',A,'  m'/ &
1515            '       vg:          ',A,'  m/s'/ &
1516            '       Gradient:    ',A,'  1/100s'/ &
1517            '       Gridpoint:   ',A)
1518425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
1519            '       Height:     ',A,'  m'/ &
1520            '       Salinity:   ',A,'  psu'/ &
1521            '       Gradient:   ',A,'  psu/100m'/ &
1522            '       Gridpoint:  ',A)
1523450 FORMAT (//' LES / Turbulence quantities:'/ &
1524              ' ---------------------------'/)
1525451 FORMAT ('   Diffusion coefficients are constant:'/ &
1526            '   Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1527452 FORMAT ('   Mixing length is limited to the Prandtl mixing lenth.')
1528453 FORMAT ('   Mixing length is limited to ',F4.2,' * z')
1529454 FORMAT ('   TKE is not allowed to fall below ',E9.2,' (m/s)**2')
1530455 FORMAT ('   initial TKE is prescribed as ',E9.2,' (m/s)**2')
1531470 FORMAT (//' Actions during the simulation:'/ &
1532              ' -----------------------------'/)
1533471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1534            '    Disturbance amplitude           :     ',F4.2, ' m/s'/       &
1535            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1536            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
1537472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1538                 ' to i/j =',I4)
1539473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1540                 1X,F5.3, ' m**2/s**2')
1541474 FORMAT ('    Random number generator used    : ',A/)
1542475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1543                 'respectively, if'/ &
1544            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1545                 ' 3D-simulation'/)
1546476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1547                 'respectively, if the'/ &
1548            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1549                 ' the 3D-simulation'/)
1550477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1551                 'respectively, if the'/ &
1552            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1553                 ' the 3D-simulation'/)
1554480 FORMAT ('    Particles:'/ &
1555            '    ---------'// &
1556            '       Particle advection is active (switched on at t = ', F7.1, &
1557                    ' s)'/ &
1558            '       Start of new particle generations every  ',F6.1,' s'/ &
1559            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
1560            '                            bottom:     ', A, ' top:         ', A/&
1561            '       Maximum particle age:                 ',F9.1,' s'/ &
1562            '       Advection stopped at t = ',F9.1,' s'/ &
1563            '       Particles are sorted every ',F9.1,' s'/)
1564481 FORMAT ('       Particles have random start positions'/)
1565482 FORMAT ('       Particles are advected only horizontally'/)
1566483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
1567484 FORMAT ('            Number of tails of the total domain: ',I10/ &
1568            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
1569            '            Maximum age of the end of the tail:  ',F8.2,' s')
1570485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
1571486 FORMAT ('       Particle statistics are written on file'/)
1572487 FORMAT ('       Number of particle groups: ',I2/)
1573488 FORMAT ('       SGS velocity components are used for particle advection'/ &
1574            '          minimum timestep for advection: ', F7.5/)
1575489 FORMAT ('       Number of particles simultaneously released at each ', &
1576                    'point: ', I5/)
1577490 FORMAT ('       Particle group ',I2,':'/ &
1578            '          Particle radius: ',E10.3, 'm')
1579491 FORMAT ('          Particle inertia is activated'/ &
1580            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
1581492 FORMAT ('          Particles are advected only passively (no inertia)'/)
1582493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
1583            '                                         y:',F8.1,' - ',F8.1,' m'/&
1584            '                                         z:',F8.1,' - ',F8.1,' m'/&
1585            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
1586                       ' m  dz = ',F8.1,' m'/)
1587494 FORMAT ('       Output of particle time series in NetCDF format every ', &
1588                    F8.2,' s'/)
1589495 FORMAT ('       Number of particles in total domain: ',I10/)
1590500 FORMAT (//' 1D-Model parameters:'/                           &
1591              ' -------------------'//                           &
1592            '    Simulation time:                   ',F8.1,' s'/ &
1593            '    Run-controll output every:         ',F8.1,' s'/ &
1594            '    Vertical profile output every:     ',F8.1,' s'/ &
1595            '    Mixing length calculation:         ',A/         &
1596            '    Dissipation calculation:           ',A/)
1597502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
1598
1599
1600 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.