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

Last change on this file since 1115 was 1115, checked in by hoffmann, 11 years ago

optimization of two-moments cloud physics

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