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

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

subdomains must have identical size, i.e. grid_matching = "match" not allowed any more
parameter grid_matching removed
some obsolete variables removed

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