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

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

further adjustments for SGI and other small changes

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