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

Last change on this file since 1127 was 1116, checked in by hoffmann, 12 years ago

last commit documented

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