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

Last change on this file since 152 was 151, checked in by raasch, 16 years ago

preliminary update for the turbulence recycling method

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