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

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

last commit documented

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