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

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

last commit documented

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