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

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

code has been put under the GNU General Public License (v3)

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