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

Last change on this file since 750 was 708, checked in by raasch, 14 years ago

last commit documented

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