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

Last change on this file since 1106 was 1106, checked in by raasch, 11 years ago

New:
---

Porting of FFT-solver for serial runs to GPU using CUDA FFT,
preprocessor lines in transpose routines rearranged, so that routines can also
be used in serial (non-parallel) mode,
transpositions also carried out in serial mode, routines fftx, fftxp replaced
by calls of fft_x, fft_x replaced by fft_x_1d in the 1D-decomposition routines
(Makefile, Makefile_check, fft_xy, poisfft, poisfft_hybrid, transpose, new: cuda_fft_interfaces)

--stdin argument for mpiexec on lckyuh, -y and -Y settings output to header (mrun)

Changed:


Module array_kind renamed precision_kind
(check_open, data_output_3d, fft_xy, modules, user_data_output_3d)

some format changes for coupled atmosphere-ocean runs (header)
small changes in code formatting (microphysics, prognostic_equations)

Errors:


bugfix: default value (0) assigned to coupling_start_time (modules)
bugfix: initial time for preruns of coupled runs is output as -coupling_start_time (data_output_profiles)

  • Property svn:keywords set to Id
File size: 76.6 KB
Line 
1 SUBROUTINE header
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2012  Leibniz University Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22! some format changes for coupled runs
23!
24! Former revisions:
25! -----------------
26! $Id: header.f90 1106 2013-03-04 05:31:38Z raasch $
27!
28! 1092 2013-02-02 11:24:22Z raasch
29! unused variables removed
30!
31! 1036 2012-10-22 13:43:42Z raasch
32! code put under GPL (PALM 3.9)
33!
34! 1031 2012-10-19 14:35:30Z raasch
35! output of netCDF data format modified
36!
37! 1015 2012-09-27 09:23:24Z raasch
38! output of Aajustment of mixing length to the Prandtl mixing length at first
39! grid point above ground removed
40!
41! 1003 2012-09-14 14:35:53Z raasch
42! output of information about equal/unequal subdomain size removed
43!
44! 1001 2012-09-13 14:08:46Z raasch
45! all actions concerning leapfrog- and upstream-spline-scheme removed
46!
47! 978 2012-08-09 08:28:32Z fricke
48! -km_damp_max, outflow_damping_width
49! +pt_damping_factor, pt_damping_width
50! +z0h
51!
52! 964 2012-07-26 09:14:24Z raasch
53! output of profil-related quantities removed
54!
55! 940 2012-07-09 14:31:00Z raasch
56! Output in case of simulations for pure neutral stratification (no pt-equation
57! solved)
58!
59! 927 2012-06-06 19:15:04Z raasch
60! output of masking_method for mg-solver
61!
62! 868 2012-03-28 12:21:07Z raasch
63! translation velocity in Galilean transformation changed to 0.6 * ug
64!
65! 833 2012-02-22 08:55:55Z maronga
66! Adjusted format for leaf area density
67!
68! 828 2012-02-21 12:00:36Z raasch
69! output of dissipation_classes + radius_classes
70!
71! 825 2012-02-19 03:03:44Z raasch
72! Output of cloud physics parameters/quantities complemented and restructured
73!
74! 767 2011-10-14 06:39:12Z raasch
75! Output of given initial u,v-profiles
76!
77! 759 2011-09-15 13:58:31Z raasch
78! output of maximum number of parallel io streams
79!
80! 707 2011-03-29 11:39:40Z raasch
81! bc_lr/ns replaced by bc_lr/ns_cyc
82!
83! 667 2010-12-23 12:06:00Z suehring/gryschka
84! Output of advection scheme.
85! Modified output of Prandtl-layer height.
86!
87! 580 2010-10-05 13:59:11Z heinze
88! Renaming of ws_vertical_gradient to subs_vertical_gradient,
89! ws_vertical_gradient_level to subs_vertical_gradient_level and
90! ws_vertical_gradient_level_ind to subs_vertical_gradient_level_i
91!
92! 493 2010-03-01 08:30:24Z raasch
93! NetCDF data output format extendend for NetCDF4/HDF5
94!
95! 449 2010-02-02 11:23:59Z raasch
96! +large scale vertical motion (subsidence/ascent)
97! Bugfix: index problem concerning gradient_level indices removed
98!
99! 410 2009-12-04 17:05:40Z letzel
100! Masked data output: + dt_domask, mask_01~20_x|y|z, mask_01~20_x|y|z_loop,
101! mask_scale|_x|y|z, masks, skip_time_domask
102!
103! 346 2009-07-06 10:13:41Z raasch
104! initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill'
105! Coupling with independent precursor runs.
106! Output of messages replaced by message handling routine.
107! Output of several additional dvr parameters
108! +canyon_height, canyon_width_x, canyon_width_y, canyon_wall_left,
109! canyon_wall_south, conserve_volume_flow_mode, dp_external, dp_level_b,
110! dp_smooth, dpdxy, u_bulk, v_bulk
111! topography_grid_convention moved from user_header
112! small bugfix concerning 3d 64bit netcdf output format
113!
114! 206 2008-10-13 14:59:11Z raasch
115! Bugfix: error in zu index in case of section_xy = -1
116!
117! 198 2008-09-17 08:55:28Z raasch
118! Format adjustments allowing output of larger revision numbers
119!
120! 197 2008-09-16 15:29:03Z raasch
121! allow 100 spectra levels instead of 10 for consistency with
122! define_netcdf_header,
123! bugfix in the output of the characteristic levels of potential temperature,
124! geostrophic wind, scalar concentration, humidity and leaf area density,
125! output of turbulence recycling informations
126!
127! 138 2007-11-28 10:03:58Z letzel
128! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
129! Allow two instead of one digit to specify isosurface and slicer variables.
130! Output of sorting frequency of particles
131!
132! 108 2007-08-24 15:10:38Z letzel
133! Output of informations for coupled model runs (boundary conditions etc.)
134! + output of momentumfluxes at the top boundary
135! Rayleigh damping for ocean, e_init
136!
137! 97 2007-06-21 08:23:15Z raasch
138! Adjustments for the ocean version.
139! use_pt_reference renamed use_reference
140!
141! 87 2007-05-22 15:46:47Z raasch
142! Bugfix: output of use_upstream_for_tke
143!
144! 82 2007-04-16 15:40:52Z raasch
145! Preprocessor strings for different linux clusters changed to "lc",
146! routine local_flush is used for buffer flushing
147!
148! 76 2007-03-29 00:58:32Z raasch
149! Output of netcdf_64bit_3d, particles-package is now part of the default code,
150! output of the loop optimization method, moisture renamed humidity,
151! output of subversion revision number
152!
153! 19 2007-02-23 04:53:48Z raasch
154! Output of scalar flux applied at top boundary
155!
156! RCS Log replace by Id keyword, revision history cleaned up
157!
158! Revision 1.63  2006/08/22 13:53:13  raasch
159! Output of dz_max
160!
161! Revision 1.1  1997/08/11 06:17:20  raasch
162! Initial revision
163!
164!
165! Description:
166! ------------
167! Writing a header with all important informations about the actual run.
168! This subroutine is called three times, two times at the beginning
169! (writing information on files RUN_CONTROL and HEADER) and one time at the
170! end of the run, then writing additional information about CPU-usage on file
171! header.
172!-----------------------------------------------------------------------------!
173
174    USE arrays_3d
175    USE control_parameters
176    USE cloud_parameters
177    USE cpulog
178    USE dvrp_variables
179    USE grid_variables
180    USE indices
181    USE model_1d
182    USE particle_attributes
183    USE pegrid
184    USE subsidence_mod
185    USE spectrum
186
187    IMPLICIT NONE
188
189    CHARACTER (LEN=1)  ::  prec
190    CHARACTER (LEN=2)  ::  do2d_mode
191    CHARACTER (LEN=5)  ::  section_chr
192    CHARACTER (LEN=10) ::  coor_chr, host_chr
193    CHARACTER (LEN=16) ::  begin_chr
194    CHARACTER (LEN=26) ::  ver_rev
195    CHARACTER (LEN=40) ::  output_format
196    CHARACTER (LEN=70) ::  char1, char2, dopr_chr, &
197                           do2d_xy, do2d_xz, do2d_yz, do3d_chr, &
198                           domask_chr, run_classification
199    CHARACTER (LEN=86) ::  coordinates, gradients, learde, slices,  &
200                           temperatures, ugcomponent, vgcomponent
201    CHARACTER (LEN=85) ::  roben, runten
202
203    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)
204
205    INTEGER ::  av, bh, blx, bly, bxl, bxr, byn, bys, ch, count, cwx, cwy,  &
206                cxl, cxr, cyn, cys, dim, i, io, j, l, ll, mpi_type
207    REAL    ::  cpuseconds_per_simulated_second
208
209!
210!-- Open the output file. At the end of the simulation, output is directed
211!-- to unit 19.
212    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
213         .NOT. simulated_time_at_begin /= simulated_time )  THEN
214       io = 15   !  header output on file RUN_CONTROL
215    ELSE
216       io = 19   !  header output on file HEADER
217    ENDIF
218    CALL check_open( io )
219
220!
221!-- At the end of the run, output file (HEADER) will be rewritten with
222!-- new informations
223    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
224
225!
226!-- Determine kind of model run
227    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
228       run_classification = '3D - restart run'
229    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
230       run_classification = '3D - run with cyclic fill of 3D - prerun data'
231    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
232       run_classification = '3D - run without 1D - prerun'
233    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
234       run_classification = '3D - run with 1D - prerun'
235    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
236       run_classification = '3D - run initialized by user'
237    ELSE
238       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
239       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
240    ENDIF
241    IF ( ocean )  THEN
242       run_classification = 'ocean - ' // run_classification
243    ELSE
244       run_classification = 'atmosphere - ' // run_classification
245    ENDIF
246
247!
248!-- Run-identification, date, time, host
249    host_chr = host(1:10)
250    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
251    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
252    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
253#if defined( __mpi2 )
254       mpi_type = 2
255#else
256       mpi_type = 1
257#endif
258       WRITE ( io, 101 )  mpi_type, coupling_mode
259    ENDIF
260    IF ( coupling_start_time /= 0.0 )  THEN
261       IF ( coupling_start_time > simulated_time_at_begin )  THEN
262          WRITE ( io, 109 )
263       ELSE
264          WRITE ( io, 114 )
265       ENDIF
266    ENDIF
267    WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr, &
268                       ADJUSTR( host_chr )
269#if defined( __parallel )
270    IF ( npex == -1  .AND.  pdims(2) /= 1 )  THEN
271       char1 = 'calculated'
272    ELSEIF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.  &
273               host(1:2) == 'lc' )  .AND.                          &
274             npex == -1  .AND.  pdims(2) == 1 )  THEN
275       char1 = 'forced'
276    ELSE
277       char1 = 'predefined'
278    ENDIF
279    IF ( threads_per_task == 1 )  THEN
280       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
281    ELSE
282       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
283                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
284    ENDIF
285    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
286           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
287         npex == -1  .AND.  pdims(2) == 1 )                      &
288    THEN
289       WRITE ( io, 106 )
290    ELSEIF ( pdims(2) == 1 )  THEN
291       WRITE ( io, 107 )  'x'
292    ELSEIF ( pdims(1) == 1 )  THEN
293       WRITE ( io, 107 )  'y'
294    ENDIF
295    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
296    IF ( numprocs /= maximum_parallel_io_streams )  THEN
297       WRITE ( io, 108 )  maximum_parallel_io_streams
298    ENDIF
299#endif
300    WRITE ( io, 99 )
301
302!
303!-- Numerical schemes
304    WRITE ( io, 110 )
305    IF ( psolver(1:7) == 'poisfft' )  THEN
306       WRITE ( io, 111 )  TRIM( fft_method )
307       IF ( psolver == 'poisfft_hybrid' )  WRITE ( io, 138 )
308    ELSEIF ( psolver == 'sor' )  THEN
309       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
310    ELSEIF ( psolver == 'multigrid' )  THEN
311       WRITE ( io, 135 )  cycle_mg, maximum_grid_level, ngsrb
312       IF ( mg_cycles == -1 )  THEN
313          WRITE ( io, 140 )  residual_limit
314       ELSE
315          WRITE ( io, 141 )  mg_cycles
316       ENDIF
317       IF ( mg_switch_to_pe0_level == 0 )  THEN
318          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
319                             nzt_mg(1)
320       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
321          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
322                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
323                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
324                             nzt_mg(mg_switch_to_pe0_level),    &
325                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
326                             nzt_mg(1)
327       ENDIF
328       IF ( masking_method )  WRITE ( io, 144 )
329    ENDIF
330    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
331    THEN
332       WRITE ( io, 142 )
333    ENDIF
334
335    IF ( momentum_advec == 'pw-scheme' )  THEN
336       WRITE ( io, 113 )
337    ELSEIF (momentum_advec == 'ws-scheme' ) THEN
338       WRITE ( io, 503 )
339    ENDIF
340    IF ( scalar_advec == 'pw-scheme' )  THEN
341       WRITE ( io, 116 )
342    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
343       WRITE ( io, 504 )
344    ELSE
345       WRITE ( io, 118 )
346    ENDIF
347
348    WRITE ( io, 139 )  TRIM( loop_optimization )
349
350    IF ( galilei_transformation )  THEN
351       IF ( use_ug_for_galilei_tr )  THEN
352          char1 = '0.6 * geostrophic wind'
353       ELSE
354          char1 = 'mean wind in model domain'
355       ENDIF
356       IF ( simulated_time_at_begin == simulated_time )  THEN
357          char2 = 'at the start of the run'
358       ELSE
359          char2 = 'at the end of the run'
360       ENDIF
361       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ), &
362                          advected_distance_x/1000.0, advected_distance_y/1000.0
363    ENDIF
364    WRITE ( io, 122 )  timestep_scheme
365    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
366    IF ( rayleigh_damping_factor /= 0.0 )  THEN
367       IF ( .NOT. ocean )  THEN
368          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
369               rayleigh_damping_factor
370       ELSE
371          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
372               rayleigh_damping_factor
373       ENDIF
374    ENDIF
375    IF ( neutral )  WRITE ( io, 131 )  pt_surface
376    IF ( humidity )  THEN
377       IF ( .NOT. cloud_physics )  THEN
378          WRITE ( io, 129 )
379       ELSE
380          WRITE ( io, 130 )
381       ENDIF
382    ENDIF
383    IF ( passive_scalar )  WRITE ( io, 134 )
384    IF ( conserve_volume_flow )  THEN
385       WRITE ( io, 150 )  conserve_volume_flow_mode
386       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
387          WRITE ( io, 151 )  u_bulk, v_bulk
388       ENDIF
389    ELSEIF ( dp_external )  THEN
390       IF ( dp_smooth )  THEN
391          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
392       ELSE
393          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
394       ENDIF
395    ENDIF
396    IF ( large_scale_subsidence )  THEN
397        WRITE ( io, 153 )
398        WRITE ( io, 154 )
399    ENDIF
400    WRITE ( io, 99 )
401
402!
403!-- Runtime and timestep informations
404    WRITE ( io, 200 )
405    IF ( .NOT. dt_fixed )  THEN
406       WRITE ( io, 201 )  dt_max, cfl_factor
407    ELSE
408       WRITE ( io, 202 )  dt
409    ENDIF
410    WRITE ( io, 203 )  simulated_time_at_begin, end_time
411
412    IF ( time_restart /= 9999999.9  .AND. &
413         simulated_time_at_begin == simulated_time )  THEN
414       IF ( dt_restart == 9999999.9 )  THEN
415          WRITE ( io, 204 )  ' Restart at:       ',time_restart
416       ELSE
417          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
418       ENDIF
419    ENDIF
420
421    IF ( simulated_time_at_begin /= simulated_time )  THEN
422       i = MAX ( log_point_s(10)%counts, 1 )
423       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0 )  THEN
424          cpuseconds_per_simulated_second = 0.0
425       ELSE
426          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
427                                            ( simulated_time -    &
428                                              simulated_time_at_begin )
429       ENDIF
430       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum, &
431                          log_point_s(10)%sum / REAL( i ),     &
432                          cpuseconds_per_simulated_second
433       IF ( time_restart /= 9999999.9  .AND.  time_restart < end_time )  THEN
434          IF ( dt_restart == 9999999.9 )  THEN
435             WRITE ( io, 204 )  ' Next restart at:     ',time_restart
436          ELSE
437             WRITE ( io, 205 )  ' Next restart at:     ',time_restart, dt_restart
438          ENDIF
439       ENDIF
440    ENDIF
441
442!
443!-- Start time for coupled runs, if independent precursor runs for atmosphere
444!-- and ocean are used or have been used. In this case, coupling_start_time
445!-- defines the time when the coupling is switched on.
446    IF ( coupling_start_time /= 0.0 )  THEN
447       WRITE ( io, 207 )  coupling_start_time
448    ENDIF
449
450!
451!-- Computational grid
452    IF ( .NOT. ocean )  THEN
453       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
454       IF ( dz_stretch_level_index < nzt+1 )  THEN
455          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
456                             dz_stretch_factor, dz_max
457       ENDIF
458    ELSE
459       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
460       IF ( dz_stretch_level_index > 0 )  THEN
461          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
462                             dz_stretch_factor, dz_max
463       ENDIF
464    ENDIF
465    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
466                       MIN( nnz+2, nzt+2 )
467    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
468
469!
470!-- Topography
471    WRITE ( io, 270 )  topography
472    SELECT CASE ( TRIM( topography ) )
473
474       CASE ( 'flat' )
475          ! no actions necessary
476
477       CASE ( 'single_building' )
478          blx = INT( building_length_x / dx )
479          bly = INT( building_length_y / dy )
480          bh  = INT( building_height / dz )
481
482          IF ( building_wall_left == 9999999.9 )  THEN
483             building_wall_left = ( nx + 1 - blx ) / 2 * dx
484          ENDIF
485          bxl = INT ( building_wall_left / dx + 0.5 )
486          bxr = bxl + blx
487
488          IF ( building_wall_south == 9999999.9 )  THEN
489             building_wall_south = ( ny + 1 - bly ) / 2 * dy
490          ENDIF
491          bys = INT ( building_wall_south / dy + 0.5 )
492          byn = bys + bly
493
494          WRITE ( io, 271 )  building_length_x, building_length_y, &
495                             building_height, bxl, bxr, bys, byn
496
497       CASE ( 'single_street_canyon' )
498          ch  = NINT( canyon_height / dz )
499          IF ( canyon_width_x /= 9999999.9 )  THEN
500!
501!--          Street canyon in y direction
502             cwx = NINT( canyon_width_x / dx )
503             IF ( canyon_wall_left == 9999999.9 )  THEN
504                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
505             ENDIF
506             cxl = NINT( canyon_wall_left / dx )
507             cxr = cxl + cwx
508             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
509
510          ELSEIF ( canyon_width_y /= 9999999.9 )  THEN
511!
512!--          Street canyon in x direction
513             cwy = NINT( canyon_width_y / dy )
514             IF ( canyon_wall_south == 9999999.9 )  THEN
515                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
516             ENDIF
517             cys = NINT( canyon_wall_south / dy )
518             cyn = cys + cwy
519             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
520          ENDIF
521
522    END SELECT
523
524    IF ( TRIM( topography ) /= 'flat' )  THEN
525       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
526          IF ( TRIM( topography ) == 'single_building' .OR.  &
527               TRIM( topography ) == 'single_street_canyon' )  THEN
528             WRITE ( io, 278 )
529          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
530             WRITE ( io, 279 )
531          ENDIF
532       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
533          WRITE ( io, 278 )
534       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
535          WRITE ( io, 279 )
536       ENDIF
537    ENDIF
538
539    IF ( plant_canopy ) THEN
540
541       WRITE ( io, 280 ) canopy_mode, pch_index, drag_coefficient
542       IF ( passive_scalar ) THEN
543          WRITE ( io, 281 ) scalar_exchange_coefficient,   &
544                            leaf_surface_concentration
545       ENDIF
546
547!
548!--    Heat flux at the top of vegetation
549       WRITE ( io, 282 ) cthf
550
551!
552!--    Leaf area density profile
553!--    Building output strings, starting with surface value
554       WRITE ( learde, '(F6.4)' )  lad_surface
555       gradients = '------'
556       slices = '     0'
557       coordinates = '   0.0'
558       i = 1
559       DO  WHILE ( lad_vertical_gradient_level_ind(i) /= -9999 )
560
561          WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
562          learde = TRIM( learde ) // ' ' // TRIM( coor_chr )
563
564          WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
565          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
566
567          WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
568          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
569
570          WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
571          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
572
573          i = i + 1
574       ENDDO
575
576       WRITE ( io, 283 )  TRIM( coordinates ), TRIM( learde ), &
577                          TRIM( gradients ), TRIM( slices )
578
579    ENDIF
580
581!
582!-- Boundary conditions
583    IF ( ibc_p_b == 0 )  THEN
584       runten = 'p(0)     = 0      |'
585    ELSEIF ( ibc_p_b == 1 )  THEN
586       runten = 'p(0)     = p(1)   |'
587    ELSE
588       runten = 'p(0)     = p(1) +R|'
589    ENDIF
590    IF ( ibc_p_t == 0 )  THEN
591       roben  = 'p(nzt+1) = 0      |'
592    ELSE
593       roben  = 'p(nzt+1) = p(nzt) |'
594    ENDIF
595
596    IF ( ibc_uv_b == 0 )  THEN
597       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
598    ELSE
599       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
600    ENDIF
601    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
602       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
603    ELSEIF ( ibc_uv_t == 0 )  THEN
604       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
605    ELSE
606       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
607    ENDIF
608
609    IF ( ibc_pt_b == 0 )  THEN
610       runten = TRIM( runten ) // ' pt(0)   = pt_surface'
611    ELSEIF ( ibc_pt_b == 1 )  THEN
612       runten = TRIM( runten ) // ' pt(0)   = pt(1)'
613    ELSEIF ( ibc_pt_b == 2 )  THEN
614       runten = TRIM( runten ) // ' pt(0) = from coupled model'
615    ENDIF
616    IF ( ibc_pt_t == 0 )  THEN
617       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
618    ELSEIF( ibc_pt_t == 1 )  THEN
619       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
620    ELSEIF( ibc_pt_t == 2 )  THEN
621       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
622
623    ENDIF
624
625    WRITE ( io, 300 )  runten, roben
626
627    IF ( .NOT. constant_diffusion )  THEN
628       IF ( ibc_e_b == 1 )  THEN
629          runten = 'e(0)     = e(1)'
630       ELSE
631          runten = 'e(0)     = e(1) = (u*/0.1)**2'
632       ENDIF
633       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
634
635       WRITE ( io, 301 )  'e', runten, roben       
636
637    ENDIF
638
639    IF ( ocean )  THEN
640       runten = 'sa(0)    = sa(1)'
641       IF ( ibc_sa_t == 0 )  THEN
642          roben =  'sa(nzt+1) = sa_surface'
643       ELSE
644          roben =  'sa(nzt+1) = sa(nzt)'
645       ENDIF
646       WRITE ( io, 301 ) 'sa', runten, roben
647    ENDIF
648
649    IF ( humidity )  THEN
650       IF ( ibc_q_b == 0 )  THEN
651          runten = 'q(0)     = q_surface'
652       ELSE
653          runten = 'q(0)     = q(1)'
654       ENDIF
655       IF ( ibc_q_t == 0 )  THEN
656          roben =  'q(nzt)   = q_top'
657       ELSE
658          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
659       ENDIF
660       WRITE ( io, 301 ) 'q', runten, roben
661    ENDIF
662
663    IF ( passive_scalar )  THEN
664       IF ( ibc_q_b == 0 )  THEN
665          runten = 's(0)     = s_surface'
666       ELSE
667          runten = 's(0)     = s(1)'
668       ENDIF
669       IF ( ibc_q_t == 0 )  THEN
670          roben =  's(nzt)   = s_top'
671       ELSE
672          roben =  's(nzt)   = s(nzt-1) + ds/dz'
673       ENDIF
674       WRITE ( io, 301 ) 's', runten, roben
675    ENDIF
676
677    IF ( use_surface_fluxes )  THEN
678       WRITE ( io, 303 )
679       IF ( constant_heatflux )  THEN
680          WRITE ( io, 306 )  surface_heatflux
681          IF ( random_heatflux )  WRITE ( io, 307 )
682       ENDIF
683       IF ( humidity  .AND.  constant_waterflux )  THEN
684          WRITE ( io, 311 ) surface_waterflux
685       ENDIF
686       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
687          WRITE ( io, 313 ) surface_waterflux
688       ENDIF
689    ENDIF
690
691    IF ( use_top_fluxes )  THEN
692       WRITE ( io, 304 )
693       IF ( coupling_mode == 'uncoupled' )  THEN
694          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
695          IF ( constant_top_heatflux )  THEN
696             WRITE ( io, 306 )  top_heatflux
697          ENDIF
698       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
699          WRITE ( io, 316 )
700       ENDIF
701       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
702          WRITE ( io, 309 )  top_salinityflux
703       ENDIF
704       IF ( humidity  .OR.  passive_scalar )  THEN
705          WRITE ( io, 315 )
706       ENDIF
707    ENDIF
708
709    IF ( prandtl_layer )  THEN
710       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length, &
711                          z0h_factor*roughness_length, kappa, &
712                          rif_min, rif_max
713       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
714       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
715          WRITE ( io, 312 )
716       ENDIF
717       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
718          WRITE ( io, 314 )
719       ENDIF
720    ELSE
721       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
722          WRITE ( io, 310 )  rif_min, rif_max
723       ENDIF
724    ENDIF
725
726    WRITE ( io, 317 )  bc_lr, bc_ns
727    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
728       WRITE ( io, 318 )  pt_damping_width, pt_damping_factor       
729       IF ( turbulent_inflow )  THEN
730          WRITE ( io, 319 )  recycling_width, recycling_plane, &
731                             inflow_damping_height, inflow_damping_width
732       ENDIF
733    ENDIF
734
735!
736!-- Listing of 1D-profiles
737    WRITE ( io, 325 )  dt_dopr_listing
738    IF ( averaging_interval_pr /= 0.0 )  THEN
739       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
740    ENDIF
741
742!
743!-- DATA output
744    WRITE ( io, 330 )
745    IF ( averaging_interval_pr /= 0.0 )  THEN
746       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
747    ENDIF
748
749!
750!-- 1D-profiles
751    dopr_chr = 'Profile:'
752    IF ( dopr_n /= 0 )  THEN
753       WRITE ( io, 331 )
754
755       output_format = ''
756       IF ( netcdf_output )  THEN
757          output_format = output_format_netcdf
758       ENDIF
759       WRITE ( io, 344 )  output_format
760
761       DO  i = 1, dopr_n
762          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
763          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
764             WRITE ( io, 332 )  dopr_chr
765             dopr_chr = '       :'
766          ENDIF
767       ENDDO
768
769       IF ( dopr_chr /= '' )  THEN
770          WRITE ( io, 332 )  dopr_chr
771       ENDIF
772       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
773       IF ( skip_time_dopr /= 0.0 )  WRITE ( io, 339 )  skip_time_dopr
774    ENDIF
775
776!
777!-- 2D-arrays
778    DO  av = 0, 1
779
780       i = 1
781       do2d_xy = ''
782       do2d_xz = ''
783       do2d_yz = ''
784       DO  WHILE ( do2d(av,i) /= ' ' )
785
786          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
787          do2d_mode = do2d(av,i)(l-1:l)
788
789          SELECT CASE ( do2d_mode )
790             CASE ( 'xy' )
791                ll = LEN_TRIM( do2d_xy )
792                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
793             CASE ( 'xz' )
794                ll = LEN_TRIM( do2d_xz )
795                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
796             CASE ( 'yz' )
797                ll = LEN_TRIM( do2d_yz )
798                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
799          END SELECT
800
801          i = i + 1
802
803       ENDDO
804
805       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
806              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
807              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) )  .AND. &
808            ( netcdf_output  .OR.  iso2d_output ) )  THEN
809
810          IF (  av == 0 )  THEN
811             WRITE ( io, 334 )  ''
812          ELSE
813             WRITE ( io, 334 )  '(time-averaged)'
814          ENDIF
815
816          IF ( do2d_at_begin )  THEN
817             begin_chr = 'and at the start'
818          ELSE
819             begin_chr = ''
820          ENDIF
821
822          output_format = ''
823          IF ( netcdf_output )  THEN
824             output_format = output_format_netcdf
825          ENDIF
826          IF ( iso2d_output )  THEN
827             IF ( netcdf_output )  THEN
828                output_format = TRIM( output_format_netcdf ) // ' and iso2d'
829             ELSE
830                output_format = 'iso2d'
831             ENDIF
832          ENDIF
833          WRITE ( io, 344 )  output_format
834
835          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
836             i = 1
837             slices = '/'
838             coordinates = '/'
839!
840!--          Building strings with index and coordinate informations of the
841!--          slices
842             DO  WHILE ( section(i,1) /= -9999 )
843
844                WRITE (section_chr,'(I5)')  section(i,1)
845                section_chr = ADJUSTL( section_chr )
846                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
847
848                IF ( section(i,1) == -1 )  THEN
849                   WRITE (coor_chr,'(F10.1)')  -1.0
850                ELSE
851                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
852                ENDIF
853                coor_chr = ADJUSTL( coor_chr )
854                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
855
856                i = i + 1
857             ENDDO
858             IF ( av == 0 )  THEN
859                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
860                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
861                                   TRIM( coordinates )
862                IF ( skip_time_do2d_xy /= 0.0 )  THEN
863                   WRITE ( io, 339 )  skip_time_do2d_xy
864                ENDIF
865             ELSE
866                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
867                                   TRIM( begin_chr ), averaging_interval, &
868                                   dt_averaging_input, 'k', TRIM( slices ), &
869                                   TRIM( coordinates )
870                IF ( skip_time_data_output_av /= 0.0 )  THEN
871                   WRITE ( io, 339 )  skip_time_data_output_av
872                ENDIF
873             ENDIF
874
875          ENDIF
876
877          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
878             i = 1
879             slices = '/'
880             coordinates = '/'
881!
882!--          Building strings with index and coordinate informations of the
883!--          slices
884             DO  WHILE ( section(i,2) /= -9999 )
885
886                WRITE (section_chr,'(I5)')  section(i,2)
887                section_chr = ADJUSTL( section_chr )
888                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
889
890                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
891                coor_chr = ADJUSTL( coor_chr )
892                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
893
894                i = i + 1
895             ENDDO
896             IF ( av == 0 )  THEN
897                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
898                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
899                                   TRIM( coordinates )
900                IF ( skip_time_do2d_xz /= 0.0 )  THEN
901                   WRITE ( io, 339 )  skip_time_do2d_xz
902                ENDIF
903             ELSE
904                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
905                                   TRIM( begin_chr ), averaging_interval, &
906                                   dt_averaging_input, 'j', TRIM( slices ), &
907                                   TRIM( coordinates )
908                IF ( skip_time_data_output_av /= 0.0 )  THEN
909                   WRITE ( io, 339 )  skip_time_data_output_av
910                ENDIF
911             ENDIF
912          ENDIF
913
914          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
915             i = 1
916             slices = '/'
917             coordinates = '/'
918!
919!--          Building strings with index and coordinate informations of the
920!--          slices
921             DO  WHILE ( section(i,3) /= -9999 )
922
923                WRITE (section_chr,'(I5)')  section(i,3)
924                section_chr = ADJUSTL( section_chr )
925                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
926
927                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
928                coor_chr = ADJUSTL( coor_chr )
929                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
930
931                i = i + 1
932             ENDDO
933             IF ( av == 0 )  THEN
934                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
935                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
936                                   TRIM( coordinates )
937                IF ( skip_time_do2d_yz /= 0.0 )  THEN
938                   WRITE ( io, 339 )  skip_time_do2d_yz
939                ENDIF
940             ELSE
941                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
942                                   TRIM( begin_chr ), averaging_interval, &
943                                   dt_averaging_input, 'i', TRIM( slices ), &
944                                   TRIM( coordinates )
945                IF ( skip_time_data_output_av /= 0.0 )  THEN
946                   WRITE ( io, 339 )  skip_time_data_output_av
947                ENDIF
948             ENDIF
949          ENDIF
950
951       ENDIF
952
953    ENDDO
954
955!
956!-- 3d-arrays
957    DO  av = 0, 1
958
959       i = 1
960       do3d_chr = ''
961       DO  WHILE ( do3d(av,i) /= ' ' )
962
963          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
964          i = i + 1
965
966       ENDDO
967
968       IF ( do3d_chr /= '' )  THEN
969          IF ( av == 0 )  THEN
970             WRITE ( io, 336 )  ''
971          ELSE
972             WRITE ( io, 336 )  '(time-averaged)'
973          ENDIF
974
975          output_format = ''
976          IF ( netcdf_output )  THEN
977             output_format = output_format_netcdf
978          ENDIF
979          IF ( avs_output )  THEN
980             IF ( netcdf_output )  THEN
981                output_format = TRIM( output_format_netcdf ) // ' and avs'
982             ELSE
983                output_format = 'avs'
984             ENDIF
985          ENDIF
986          WRITE ( io, 344 )  output_format
987
988          IF ( do3d_at_begin )  THEN
989             begin_chr = 'and at the start'
990          ELSE
991             begin_chr = ''
992          ENDIF
993          IF ( av == 0 )  THEN
994             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
995                                zu(nz_do3d), nz_do3d
996          ELSE
997             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
998                                TRIM( begin_chr ), averaging_interval, &
999                                dt_averaging_input, zu(nz_do3d), nz_do3d
1000          ENDIF
1001
1002          IF ( do3d_compress )  THEN
1003             do3d_chr = ''
1004             i = 1
1005             DO WHILE ( do3d(av,i) /= ' ' )
1006
1007                SELECT CASE ( do3d(av,i) )
1008                   CASE ( 'u' )
1009                      j = 1
1010                   CASE ( 'v' )
1011                      j = 2
1012                   CASE ( 'w' )
1013                      j = 3
1014                   CASE ( 'p' )
1015                      j = 4
1016                   CASE ( 'pt' )
1017                      j = 5
1018                END SELECT
1019                WRITE ( prec, '(I1)' )  plot_3d_precision(j)%precision
1020                do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // &
1021                           ':' // prec // ','
1022                i = i + 1
1023
1024             ENDDO
1025             WRITE ( io, 338 )  do3d_chr
1026
1027          ENDIF
1028
1029          IF ( av == 0 )  THEN
1030             IF ( skip_time_do3d /= 0.0 )  THEN
1031                WRITE ( io, 339 )  skip_time_do3d
1032             ENDIF
1033          ELSE
1034             IF ( skip_time_data_output_av /= 0.0 )  THEN
1035                WRITE ( io, 339 )  skip_time_data_output_av
1036             ENDIF
1037          ENDIF
1038
1039       ENDIF
1040
1041    ENDDO
1042
1043!
1044!-- masked arrays
1045    IF ( masks > 0 )  WRITE ( io, 345 )  &
1046         mask_scale_x, mask_scale_y, mask_scale_z
1047    DO  mid = 1, masks
1048       DO  av = 0, 1
1049
1050          i = 1
1051          domask_chr = ''
1052          DO  WHILE ( domask(mid,av,i) /= ' ' )
1053             domask_chr = TRIM( domask_chr ) // ' ' //  &
1054                          TRIM( domask(mid,av,i) ) // ','
1055             i = i + 1
1056          ENDDO
1057
1058          IF ( domask_chr /= '' )  THEN
1059             IF ( av == 0 )  THEN
1060                WRITE ( io, 346 )  '', mid
1061             ELSE
1062                WRITE ( io, 346 )  ' (time-averaged)', mid
1063             ENDIF
1064
1065             output_format = ' '
1066             IF ( netcdf_output )  THEN
1067                output_format = output_format_netcdf
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          output_format = output_format_netcdf
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          output_format = output_format_netcdf
1188       ENDIF
1189       WRITE ( io, 344 )  output_format
1190       WRITE ( io, 371 )  dt_dosp
1191       IF ( skip_time_dosp /= 0.0 )  WRITE ( io, 339 )  skip_time_dosp
1192       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
1193                          ( spectra_direction(i), i = 1,10 ),  &
1194                          ( comp_spectra_level(i), i = 1,100 ), &
1195                          ( plot_spectra_level(i), i = 1,100 ), &
1196                          averaging_interval_sp, dt_averaging_input_pr
1197    ENDIF
1198#endif
1199
1200    WRITE ( io, 99 )
1201
1202!
1203!-- Physical quantities
1204    WRITE ( io, 400 )
1205
1206!
1207!-- Geostrophic parameters
1208    WRITE ( io, 410 )  omega, phi, f, fs
1209
1210!
1211!-- Other quantities
1212    WRITE ( io, 411 )  g
1213    IF ( use_reference )  THEN
1214       IF ( ocean )  THEN
1215          WRITE ( io, 412 )  prho_reference
1216       ELSE
1217          WRITE ( io, 413 )  pt_reference
1218       ENDIF
1219    ENDIF
1220
1221!
1222!-- Cloud physics parameters
1223    IF ( cloud_physics ) THEN
1224       WRITE ( io, 415 )
1225       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
1226    ENDIF
1227
1228!-- Profile of the geostrophic wind (component ug)
1229!-- Building output strings
1230    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
1231    gradients = '------'
1232    slices = '     0'
1233    coordinates = '   0.0'
1234    i = 1
1235    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
1236     
1237       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
1238       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
1239
1240       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
1241       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1242
1243       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
1244       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1245
1246       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
1247       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1248
1249       IF ( i == 10 )  THEN
1250          EXIT
1251       ELSE
1252          i = i + 1
1253       ENDIF
1254
1255    ENDDO
1256
1257    WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
1258                       TRIM( gradients ), TRIM( slices )
1259
1260!-- Profile of the geostrophic wind (component vg)
1261!-- Building output strings
1262    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
1263    gradients = '------'
1264    slices = '     0'
1265    coordinates = '   0.0'
1266    i = 1
1267    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
1268
1269       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
1270       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
1271
1272       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
1273       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1274
1275       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
1276       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1277
1278       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
1279       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1280
1281       IF ( i == 10 )  THEN
1282          EXIT
1283       ELSE
1284          i = i + 1
1285       ENDIF
1286 
1287    ENDDO
1288
1289    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
1290                       TRIM( gradients ), TRIM( slices )
1291
1292!
1293!-- Initial wind profiles
1294    IF ( u_profile(1) /= 9999999.9 )  WRITE ( io, 427 )
1295
1296!
1297!-- Initial temperature profile
1298!-- Building output strings, starting with surface temperature
1299    WRITE ( temperatures, '(F6.2)' )  pt_surface
1300    gradients = '------'
1301    slices = '     0'
1302    coordinates = '   0.0'
1303    i = 1
1304    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1305
1306       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1307       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1308
1309       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1310       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1311
1312       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1313       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1314
1315       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1316       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1317
1318       IF ( i == 10 )  THEN
1319          EXIT
1320       ELSE
1321          i = i + 1
1322       ENDIF
1323
1324    ENDDO
1325
1326    WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1327                       TRIM( gradients ), TRIM( slices )
1328
1329!
1330!-- Initial humidity profile
1331!-- Building output strings, starting with surface humidity
1332    IF ( humidity  .OR.  passive_scalar )  THEN
1333       WRITE ( temperatures, '(E8.1)' )  q_surface
1334       gradients = '--------'
1335       slices = '       0'
1336       coordinates = '     0.0'
1337       i = 1
1338       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1339         
1340          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1341          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1342
1343          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1344          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1345         
1346          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1347          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1348         
1349          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1350          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1351
1352          IF ( i == 10 )  THEN
1353             EXIT
1354          ELSE
1355             i = i + 1
1356          ENDIF
1357
1358       ENDDO
1359
1360       IF ( humidity )  THEN
1361          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1362                             TRIM( gradients ), TRIM( slices )
1363       ELSE
1364          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1365                             TRIM( gradients ), TRIM( slices )
1366       ENDIF
1367    ENDIF
1368
1369!
1370!-- Initial salinity profile
1371!-- Building output strings, starting with surface salinity
1372    IF ( ocean )  THEN
1373       WRITE ( temperatures, '(F6.2)' )  sa_surface
1374       gradients = '------'
1375       slices = '     0'
1376       coordinates = '   0.0'
1377       i = 1
1378       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1379
1380          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1381          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1382
1383          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1384          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1385
1386          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1387          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1388
1389          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1390          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1391
1392          IF ( i == 10 )  THEN
1393             EXIT
1394          ELSE
1395             i = i + 1
1396          ENDIF
1397
1398       ENDDO
1399
1400       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1401                          TRIM( gradients ), TRIM( slices )
1402    ENDIF
1403
1404!
1405!-- Profile for the large scale vertial velocity
1406!-- Building output strings, starting with surface value
1407    IF ( large_scale_subsidence )  THEN
1408       temperatures = '   0.0'
1409       gradients = '------'
1410       slices = '     0'
1411       coordinates = '   0.0'
1412       i = 1
1413       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
1414
1415          WRITE (coor_chr,'(E10.2,7X)')  &
1416                                w_subs(subs_vertical_gradient_level_i(i))
1417          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1418
1419          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
1420          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1421
1422          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
1423          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1424
1425          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
1426          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1427
1428          IF ( i == 10 )  THEN
1429             EXIT
1430          ELSE
1431             i = i + 1
1432          ENDIF
1433
1434       ENDDO
1435
1436       WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
1437                          TRIM( gradients ), TRIM( slices )
1438    ENDIF
1439
1440!
1441!-- Cloud physcis parameters / quantities / numerical methods
1442    WRITE ( io, 430 )
1443    IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets)  THEN
1444       WRITE ( io, 431 )
1445    ELSEIF ( humidity  .AND.  cloud_physics )  THEN
1446       WRITE ( io, 432 )
1447       IF ( radiation )      WRITE ( io, 132 )
1448       IF ( precipitation )  WRITE ( io, 133 )
1449    ELSEIF ( humidity  .AND.  cloud_droplets )  THEN
1450       WRITE ( io, 433 )
1451       IF ( curvature_solution_effects )  WRITE ( io, 434 )
1452       IF ( collision_kernel /= 'none' )  THEN
1453          WRITE ( io, 435 )  TRIM( collision_kernel )
1454          IF ( collision_kernel(6:9) == 'fast' )  THEN
1455             WRITE ( io, 436 )  radius_classes, dissipation_classes
1456          ENDIF
1457       ELSE
1458          WRITE ( io, 437 )
1459       ENDIF
1460    ENDIF
1461
1462!
1463!-- LES / turbulence parameters
1464    WRITE ( io, 450 )
1465
1466!--
1467! ... LES-constants used must still be added here
1468!--
1469    IF ( constant_diffusion )  THEN
1470       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1471                          prandtl_number
1472    ENDIF
1473    IF ( .NOT. constant_diffusion)  THEN
1474       IF ( e_init > 0.0 )  WRITE ( io, 455 )  e_init
1475       IF ( e_min > 0.0 )  WRITE ( io, 454 )  e_min
1476       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1477    ENDIF
1478
1479!
1480!-- Special actions during the run
1481    WRITE ( io, 470 )
1482    IF ( create_disturbances )  THEN
1483       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1484                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1485                          zu(disturbance_level_ind_t), disturbance_level_ind_t
1486       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
1487          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1488       ELSE
1489          WRITE ( io, 473 )  disturbance_energy_limit
1490       ENDIF
1491       WRITE ( io, 474 )  TRIM( random_generator )
1492    ENDIF
1493    IF ( pt_surface_initial_change /= 0.0 )  THEN
1494       WRITE ( io, 475 )  pt_surface_initial_change
1495    ENDIF
1496    IF ( humidity  .AND.  q_surface_initial_change /= 0.0 )  THEN
1497       WRITE ( io, 476 )  q_surface_initial_change       
1498    ENDIF
1499    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0 )  THEN
1500       WRITE ( io, 477 )  q_surface_initial_change       
1501    ENDIF
1502
1503    IF ( particle_advection )  THEN
1504!
1505!--    Particle attributes
1506       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1507                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
1508                          end_time_prel, dt_sort_particles
1509       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1510       IF ( random_start_position )  WRITE ( io, 481 )
1511       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1512       WRITE ( io, 495 )  total_number_of_particles
1513       IF ( use_particle_tails  .AND.  maximum_number_of_tailpoints /= 0 )  THEN
1514          WRITE ( io, 483 )  maximum_number_of_tailpoints
1515          IF ( minimum_tailpoint_distance /= 0 )  THEN
1516             WRITE ( io, 484 )  total_number_of_tails,      &
1517                                minimum_tailpoint_distance, &
1518                                maximum_tailpoint_age
1519          ENDIF
1520       ENDIF
1521       IF ( dt_write_particle_data /= 9999999.9 )  THEN
1522          WRITE ( io, 485 )  dt_write_particle_data
1523          output_format = ' '
1524          IF ( netcdf_output )  THEN
1525             IF ( netcdf_data_format > 1 )  THEN
1526                output_format = 'netcdf (64 bit offset) and binary'
1527             ELSE
1528                output_format = 'netcdf and binary'
1529             ENDIF
1530          ELSE
1531             output_format = 'binary'
1532          ENDIF
1533          WRITE ( io, 344 )  output_format
1534       ENDIF
1535       IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1536       IF ( write_particle_statistics )  WRITE ( io, 486 )
1537
1538       WRITE ( io, 487 )  number_of_particle_groups
1539
1540       DO  i = 1, number_of_particle_groups
1541          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9 )  THEN
1542             WRITE ( io, 490 )  i, 0.0
1543             WRITE ( io, 492 )
1544          ELSE
1545             WRITE ( io, 490 )  i, radius(i)
1546             IF ( density_ratio(i) /= 0.0 )  THEN
1547                WRITE ( io, 491 )  density_ratio(i)
1548             ELSE
1549                WRITE ( io, 492 )
1550             ENDIF
1551          ENDIF
1552          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1553                             pdx(i), pdy(i), pdz(i)
1554          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
1555       ENDDO
1556
1557    ENDIF
1558
1559
1560!
1561!-- Parameters of 1D-model
1562    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1563       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1564                          mixing_length_1d, dissipation_1d
1565       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1566          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1567       ENDIF
1568    ENDIF
1569
1570!
1571!-- User-defined informations
1572    CALL user_header( io )
1573
1574    WRITE ( io, 99 )
1575
1576!
1577!-- Write buffer contents to disc immediately
1578    CALL local_flush( io )
1579
1580!
1581!-- Here the FORMATs start
1582
1583 99 FORMAT (1X,78('-'))
1584100 FORMAT (/1X,'******************************',6X,42('-')/        &
1585            1X,'* ',A,' *',6X,A/                               &
1586            1X,'******************************',6X,42('-'))
1587101 FORMAT (37X,'coupled run using MPI-',I1,': ',A/ &
1588            37X,42('-'))
1589102 FORMAT (/' Date:                 ',A8,6X,'Run:       ',A20/      &
1590            ' Time:                 ',A8,6X,'Run-No.:   ',I2.2/     &
1591            ' Run on host:        ',A10)
1592#if defined( __parallel )
1593103 FORMAT (' Number of PEs:',10X,I6,6X,'Processor grid (x,y): (',I3,',',I3, &
1594              ')',1X,A)
1595104 FORMAT (' Number of PEs:',8X,I5,9X,'Tasks:',I4,'   threads per task:',I4/ &
1596              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
1597105 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1598106 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
1599            37X,'because the job is running on an SMP-cluster')
1600107 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
1601108 FORMAT (37X,'Max. # of parallel I/O streams is ',I5)
1602109 FORMAT (37X,'Precursor run for coupled atmos-ocean run'/ &
1603            37X,42('-'))
1604114 FORMAT (37X,'Coupled atmosphere-ocean run following'/ &
1605            37X,'independent precursor runs'/             &
1606            37X,42('-'))
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 ( ' 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)
1935453 FORMAT ('    Mixing length is limited to ',F4.2,' * z')
1936454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
1937455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
1938470 FORMAT (//' Actions during the simulation:'/ &
1939              ' -----------------------------'/)
1940471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1941            '    Disturbance amplitude           :     ',F4.2, ' m/s'/       &
1942            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1943            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
1944472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1945                 ' to i/j =',I4)
1946473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1947                 1X,F5.3, ' m**2/s**2')
1948474 FORMAT ('    Random number generator used    : ',A/)
1949475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1950                 'respectively, if'/ &
1951            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1952                 ' 3D-simulation'/)
1953476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1954                 'respectively, if the'/ &
1955            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1956                 ' the 3D-simulation'/)
1957477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1958                 'respectively, if the'/ &
1959            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1960                 ' the 3D-simulation'/)
1961480 FORMAT ('    Particles:'/ &
1962            '    ---------'// &
1963            '       Particle advection is active (switched on at t = ', F7.1, &
1964                    ' s)'/ &
1965            '       Start of new particle generations every  ',F6.1,' s'/ &
1966            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
1967            '                            bottom:     ', A, ' top:         ', A/&
1968            '       Maximum particle age:                 ',F9.1,' s'/ &
1969            '       Advection stopped at t = ',F9.1,' s'/ &
1970            '       Particles are sorted every ',F9.1,' s'/)
1971481 FORMAT ('       Particles have random start positions'/)
1972482 FORMAT ('          Particles are advected only horizontally'/)
1973483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
1974484 FORMAT ('            Number of tails of the total domain: ',I10/ &
1975            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
1976            '            Maximum age of the end of the tail:  ',F8.2,' s')
1977485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
1978486 FORMAT ('       Particle statistics are written on file'/)
1979487 FORMAT ('       Number of particle groups: ',I2/)
1980488 FORMAT ('       SGS velocity components are used for particle advection'/ &
1981            '          minimum timestep for advection: ', F7.5/)
1982489 FORMAT ('       Number of particles simultaneously released at each ', &
1983                    'point: ', I5/)
1984490 FORMAT ('       Particle group ',I2,':'/ &
1985            '          Particle radius: ',E10.3, 'm')
1986491 FORMAT ('          Particle inertia is activated'/ &
1987            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
1988492 FORMAT ('          Particles are advected only passively (no inertia)'/)
1989493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
1990            '                                         y:',F8.1,' - ',F8.1,' m'/&
1991            '                                         z:',F8.1,' - ',F8.1,' m'/&
1992            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
1993                       ' m  dz = ',F8.1,' m'/)
1994494 FORMAT ('       Output of particle time series in NetCDF format every ', &
1995                    F8.2,' s'/)
1996495 FORMAT ('       Number of particles in total domain: ',I10/)
1997500 FORMAT (//' 1D-Model parameters:'/                           &
1998              ' -------------------'//                           &
1999            '    Simulation time:                   ',F8.1,' s'/ &
2000            '    Run-controll output every:         ',F8.1,' s'/ &
2001            '    Vertical profile output every:     ',F8.1,' s'/ &
2002            '    Mixing length calculation:         ',A/         &
2003            '    Dissipation calculation:           ',A/)
2004502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
2005503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2006504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
2007
2008
2009 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.