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

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

netCDF4 without parallel file support implemented
additional define string netcdf4_parallel is required to switch on parallel file support

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