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

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

last commit documented

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