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

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

tridia-solver moved to seperate module; the tridiagonal matrix coefficients of array tri are calculated only once at the beginning

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