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

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

file headers updated for the next release 3.5

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