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

Last change on this file since 829 was 829, checked in by raasch, 12 years ago

last commit documented

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