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

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

last commit documented

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