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

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

Starting with changes required for GPU optimization. OpenACC statements for using NVIDIA GPUs added.
Adjustment of mixing length to the Prandtl mixing length at first grid point above ground removed.
mask array is set zero for ghost boundaries

  • Property svn:keywords set to Id
File size: 77.0 KB
Line 
1 SUBROUTINE header
2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
6! output of Aajustment of mixing length to the Prandtl mixing length at first
7! grid point above ground removed
8!
9! Former revisions:
10! -----------------
11! $Id: header.f90 1015 2012-09-27 09:23:24Z raasch $
12!
13! 1003 2012-09-14 14:35:53Z raasch
14! output of information about equal/unequal subdomain size removed
15!
16! 1001 2012-09-13 14:08:46Z raasch
17! all actions concerning leapfrog- and upstream-spline-scheme removed
18!
19! 978 2012-08-09 08:28:32Z fricke
20! -km_damp_max, outflow_damping_width
21! +pt_damping_factor, pt_damping_width
22! +z0h
23!
24! 964 2012-07-26 09:14:24Z raasch
25! output of profil-related quantities removed
26!
27! 940 2012-07-09 14:31:00Z raasch
28! Output in case of simulations for pure neutral stratification (no pt-equation
29! solved)
30!
31! 927 2012-06-06 19:15:04Z raasch
32! output of masking_method for mg-solver
33!
34! 868 2012-03-28 12:21:07Z raasch
35! translation velocity in Galilean transformation changed to 0.6 * ug
36!
37! 833 2012-02-22 08:55:55Z maronga
38! Adjusted format for leaf area density
39!
40! 828 2012-02-21 12:00:36Z raasch
41! output of dissipation_classes + radius_classes
42!
43! 825 2012-02-19 03:03:44Z raasch
44! Output of cloud physics parameters/quantities complemented and restructured
45!
46! 767 2011-10-14 06:39:12Z raasch
47! Output of given initial u,v-profiles
48!
49! 759 2011-09-15 13:58:31Z raasch
50! output of maximum number of parallel io streams
51!
52! 707 2011-03-29 11:39:40Z raasch
53! bc_lr/ns replaced by bc_lr/ns_cyc
54!
55! 667 2010-12-23 12:06:00Z suehring/gryschka
56! Output of advection scheme.
57! Modified output of Prandtl-layer height.
58!
59! 580 2010-10-05 13:59:11Z heinze
60! Renaming of ws_vertical_gradient to subs_vertical_gradient,
61! ws_vertical_gradient_level to subs_vertical_gradient_level and
62! ws_vertical_gradient_level_ind to subs_vertical_gradient_level_i
63!
64! 493 2010-03-01 08:30:24Z raasch
65! NetCDF data output format extendend for NetCDF4/HDF5
66!
67! 449 2010-02-02 11:23:59Z raasch
68! +large scale vertical motion (subsidence/ascent)
69! Bugfix: index problem concerning gradient_level indices removed
70!
71! 410 2009-12-04 17:05:40Z letzel
72! Masked data output: + dt_domask, mask_01~20_x|y|z, mask_01~20_x|y|z_loop,
73! mask_scale|_x|y|z, masks, skip_time_domask
74!
75! 346 2009-07-06 10:13:41Z raasch
76! initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill'
77! Coupling with independent precursor runs.
78! Output of messages replaced by message handling routine.
79! Output of several additional dvr parameters
80! +canyon_height, canyon_width_x, canyon_width_y, canyon_wall_left,
81! canyon_wall_south, conserve_volume_flow_mode, dp_external, dp_level_b,
82! dp_smooth, dpdxy, u_bulk, v_bulk
83! topography_grid_convention moved from user_header
84! small bugfix concerning 3d 64bit netcdf output format
85!
86! 206 2008-10-13 14:59:11Z raasch
87! Bugfix: error in zu index in case of section_xy = -1
88!
89! 198 2008-09-17 08:55:28Z raasch
90! Format adjustments allowing output of larger revision numbers
91!
92! 197 2008-09-16 15:29:03Z raasch
93! allow 100 spectra levels instead of 10 for consistency with
94! define_netcdf_header,
95! bugfix in the output of the characteristic levels of potential temperature,
96! geostrophic wind, scalar concentration, humidity and leaf area density,
97! output of turbulence recycling informations
98!
99! 138 2007-11-28 10:03:58Z letzel
100! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
101! Allow two instead of one digit to specify isosurface and slicer variables.
102! Output of sorting frequency of particles
103!
104! 108 2007-08-24 15:10:38Z letzel
105! Output of informations for coupled model runs (boundary conditions etc.)
106! + output of momentumfluxes at the top boundary
107! Rayleigh damping for ocean, e_init
108!
109! 97 2007-06-21 08:23:15Z raasch
110! Adjustments for the ocean version.
111! use_pt_reference renamed use_reference
112!
113! 87 2007-05-22 15:46:47Z raasch
114! Bugfix: output of use_upstream_for_tke
115!
116! 82 2007-04-16 15:40:52Z raasch
117! Preprocessor strings for different linux clusters changed to "lc",
118! routine local_flush is used for buffer flushing
119!
120! 76 2007-03-29 00:58:32Z raasch
121! Output of netcdf_64bit_3d, particles-package is now part of the default code,
122! output of the loop optimization method, moisture renamed humidity,
123! output of subversion revision number
124!
125! 19 2007-02-23 04:53:48Z raasch
126! Output of scalar flux applied at top boundary
127!
128! RCS Log replace by Id keyword, revision history cleaned up
129!
130! Revision 1.63  2006/08/22 13:53:13  raasch
131! Output of dz_max
132!
133! Revision 1.1  1997/08/11 06:17:20  raasch
134! Initial revision
135!
136!
137! Description:
138! ------------
139! Writing a header with all important informations about the actual run.
140! This subroutine is called three times, two times at the beginning
141! (writing information on files RUN_CONTROL and HEADER) and one time at the
142! end of the run, then writing additional information about CPU-usage on file
143! header.
144!-----------------------------------------------------------------------------!
145
146    USE arrays_3d
147    USE control_parameters
148    USE cloud_parameters
149    USE cpulog
150    USE dvrp_variables
151    USE grid_variables
152    USE indices
153    USE model_1d
154    USE particle_attributes
155    USE pegrid
156    USE subsidence_mod
157    USE spectrum
158
159    IMPLICIT NONE
160
161    CHARACTER (LEN=1)  ::  prec
162    CHARACTER (LEN=2)  ::  do2d_mode
163    CHARACTER (LEN=5)  ::  section_chr
164    CHARACTER (LEN=9)  ::  time_to_string
165    CHARACTER (LEN=10) ::  coor_chr, host_chr
166    CHARACTER (LEN=16) ::  begin_chr
167    CHARACTER (LEN=23) ::  ver_rev
168    CHARACTER (LEN=40) ::  output_format
169    CHARACTER (LEN=70) ::  char1, char2, dopr_chr, &
170                           do2d_xy, do2d_xz, do2d_yz, do3d_chr, &
171                           domask_chr, run_classification
172    CHARACTER (LEN=86) ::  coordinates, gradients, learde, slices,  &
173                           temperatures, ugcomponent, vgcomponent
174    CHARACTER (LEN=85) ::  roben, runten
175
176    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)
177
178    INTEGER ::  av, bh, blx, bly, bxl, bxr, byn, bys, ch, count, cwx, cwy,  &
179         cxl, cxr, cyn, cys, dim, i, ihost, io, j, l, ll, m, mpi_type
180    REAL    ::  cpuseconds_per_simulated_second
181
182!
183!-- Open the output file. At the end of the simulation, output is directed
184!-- to unit 19.
185    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
186         .NOT. simulated_time_at_begin /= simulated_time )  THEN
187       io = 15   !  header output on file RUN_CONTROL
188    ELSE
189       io = 19   !  header output on file HEADER
190    ENDIF
191    CALL check_open( io )
192
193!
194!-- At the end of the run, output file (HEADER) will be rewritten with
195!-- new informations
196    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
197
198!
199!-- Determine kind of model run
200    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
201       run_classification = '3D - restart run'
202    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
203       run_classification = '3D - run with cyclic fill of 3D - prerun data'
204    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
205       run_classification = '3D - run without 1D - prerun'
206    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
207       run_classification = '3D - run with 1D - prerun'
208    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
209       run_classification = '3D - run initialized by user'
210    ELSE
211       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
212       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
213    ENDIF
214    IF ( ocean )  THEN
215       run_classification = 'ocean - ' // run_classification
216    ELSE
217       run_classification = 'atmosphere - ' // run_classification
218    ENDIF
219
220!
221!-- Run-identification, date, time, host
222    host_chr = host(1:10)
223    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
224    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
225    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
226#if defined( __mpi2 )
227       mpi_type = 2
228#else
229       mpi_type = 1
230#endif
231       WRITE ( io, 101 )  mpi_type, coupling_mode
232    ENDIF
233    WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr, &
234                       ADJUSTR( host_chr )
235#if defined( __parallel )
236    IF ( npex == -1  .AND.  pdims(2) /= 1 )  THEN
237       char1 = 'calculated'
238    ELSEIF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.  &
239               host(1:2) == 'lc' )  .AND.                          &
240             npex == -1  .AND.  pdims(2) == 1 )  THEN
241       char1 = 'forced'
242    ELSE
243       char1 = 'predefined'
244    ENDIF
245    IF ( threads_per_task == 1 )  THEN
246       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
247    ELSE
248       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
249                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
250    ENDIF
251    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
252           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
253         npex == -1  .AND.  pdims(2) == 1 )                      &
254    THEN
255       WRITE ( io, 106 )
256    ELSEIF ( pdims(2) == 1 )  THEN
257       WRITE ( io, 107 )  'x'
258    ELSEIF ( pdims(1) == 1 )  THEN
259       WRITE ( io, 107 )  'y'
260    ENDIF
261    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
262    IF ( numprocs /= maximum_parallel_io_streams )  THEN
263       WRITE ( io, 108 )  maximum_parallel_io_streams
264    ENDIF
265#endif
266    WRITE ( io, 99 )
267
268!
269!-- Numerical schemes
270    WRITE ( io, 110 )
271    IF ( psolver(1:7) == 'poisfft' )  THEN
272       WRITE ( io, 111 )  TRIM( fft_method )
273       IF ( psolver == 'poisfft_hybrid' )  WRITE ( io, 138 )
274    ELSEIF ( psolver == 'sor' )  THEN
275       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
276    ELSEIF ( psolver == 'multigrid' )  THEN
277       WRITE ( io, 135 )  cycle_mg, maximum_grid_level, ngsrb
278       IF ( mg_cycles == -1 )  THEN
279          WRITE ( io, 140 )  residual_limit
280       ELSE
281          WRITE ( io, 141 )  mg_cycles
282       ENDIF
283       IF ( mg_switch_to_pe0_level == 0 )  THEN
284          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
285                             nzt_mg(1)
286       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
287          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
288                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
289                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
290                             nzt_mg(mg_switch_to_pe0_level),    &
291                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
292                             nzt_mg(1)
293       ENDIF
294       IF ( masking_method )  WRITE ( io, 144 )
295    ENDIF
296    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
297    THEN
298       WRITE ( io, 142 )
299    ENDIF
300
301    IF ( momentum_advec == 'pw-scheme' )  THEN
302       WRITE ( io, 113 )
303    ELSEIF (momentum_advec == 'ws-scheme' ) THEN
304       WRITE ( io, 503 )
305    ENDIF
306    IF ( scalar_advec == 'pw-scheme' )  THEN
307       WRITE ( io, 116 )
308    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
309       WRITE ( io, 504 )
310    ELSE
311       WRITE ( io, 118 )
312    ENDIF
313
314    WRITE ( io, 139 )  TRIM( loop_optimization )
315
316    IF ( galilei_transformation )  THEN
317       IF ( use_ug_for_galilei_tr )  THEN
318          char1 = '0.6 * geostrophic wind'
319       ELSE
320          char1 = 'mean wind in model domain'
321       ENDIF
322       IF ( simulated_time_at_begin == simulated_time )  THEN
323          char2 = 'at the start of the run'
324       ELSE
325          char2 = 'at the end of the run'
326       ENDIF
327       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ), &
328                          advected_distance_x/1000.0, advected_distance_y/1000.0
329    ENDIF
330    WRITE ( io, 122 )  timestep_scheme
331    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
332    IF ( rayleigh_damping_factor /= 0.0 )  THEN
333       IF ( .NOT. ocean )  THEN
334          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
335               rayleigh_damping_factor
336       ELSE
337          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
338               rayleigh_damping_factor
339       ENDIF
340    ENDIF
341    IF ( neutral )  WRITE ( io, 131 )  pt_surface
342    IF ( humidity )  THEN
343       IF ( .NOT. cloud_physics )  THEN
344          WRITE ( io, 129 )
345       ELSE
346          WRITE ( io, 130 )
347       ENDIF
348    ENDIF
349    IF ( passive_scalar )  WRITE ( io, 134 )
350    IF ( conserve_volume_flow )  THEN
351       WRITE ( io, 150 )  conserve_volume_flow_mode
352       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
353          WRITE ( io, 151 )  u_bulk, v_bulk
354       ENDIF
355    ELSEIF ( dp_external )  THEN
356       IF ( dp_smooth )  THEN
357          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
358       ELSE
359          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
360       ENDIF
361    ENDIF
362    IF ( large_scale_subsidence )  THEN
363        WRITE ( io, 153 )
364        WRITE ( io, 154 )
365    ENDIF
366    WRITE ( io, 99 )
367
368!
369!-- Runtime and timestep informations
370    WRITE ( io, 200 )
371    IF ( .NOT. dt_fixed )  THEN
372       WRITE ( io, 201 )  dt_max, cfl_factor
373    ELSE
374       WRITE ( io, 202 )  dt
375    ENDIF
376    WRITE ( io, 203 )  simulated_time_at_begin, end_time
377
378    IF ( time_restart /= 9999999.9  .AND. &
379         simulated_time_at_begin == simulated_time )  THEN
380       IF ( dt_restart == 9999999.9 )  THEN
381          WRITE ( io, 204 )  ' Restart at:       ',time_restart
382       ELSE
383          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
384       ENDIF
385    ENDIF
386
387    IF ( simulated_time_at_begin /= simulated_time )  THEN
388       i = MAX ( log_point_s(10)%counts, 1 )
389       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0 )  THEN
390          cpuseconds_per_simulated_second = 0.0
391       ELSE
392          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
393                                            ( simulated_time -    &
394                                              simulated_time_at_begin )
395       ENDIF
396       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum, &
397                          log_point_s(10)%sum / REAL( i ),     &
398                          cpuseconds_per_simulated_second
399       IF ( time_restart /= 9999999.9  .AND.  time_restart < end_time )  THEN
400          IF ( dt_restart == 9999999.9 )  THEN
401             WRITE ( io, 204 )  ' Next restart at:  ',time_restart
402          ELSE
403             WRITE ( io, 205 )  ' Next restart at:  ',time_restart, dt_restart
404          ENDIF
405       ENDIF
406    ENDIF
407
408!
409!-- Start time for coupled runs, if independent precursor runs for atmosphere
410!-- and ocean are used. In this case, coupling_start_time defines the time
411!-- when the coupling is switched on.
412    IF ( coupling_start_time /= 0.0 )  THEN
413       IF ( coupling_start_time >= simulated_time_at_begin )  THEN
414          char1 = 'Precursor run for a coupled atmosphere-ocean run'
415       ELSE
416          char1 = 'Coupled atmosphere-ocean run following independent ' // &
417                  'precursor runs'
418       ENDIF
419       WRITE ( io, 207 )  char1, coupling_start_time
420    ENDIF
421
422!
423!-- Computational grid
424    IF ( .NOT. ocean )  THEN
425       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
426       IF ( dz_stretch_level_index < nzt+1 )  THEN
427          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
428                             dz_stretch_factor, dz_max
429       ENDIF
430    ELSE
431       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
432       IF ( dz_stretch_level_index > 0 )  THEN
433          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
434                             dz_stretch_factor, dz_max
435       ENDIF
436    ENDIF
437    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
438                       MIN( nnz+2, nzt+2 )
439    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
440
441!
442!-- Topography
443    WRITE ( io, 270 )  topography
444    SELECT CASE ( TRIM( topography ) )
445
446       CASE ( 'flat' )
447          ! no actions necessary
448
449       CASE ( 'single_building' )
450          blx = INT( building_length_x / dx )
451          bly = INT( building_length_y / dy )
452          bh  = INT( building_height / dz )
453
454          IF ( building_wall_left == 9999999.9 )  THEN
455             building_wall_left = ( nx + 1 - blx ) / 2 * dx
456          ENDIF
457          bxl = INT ( building_wall_left / dx + 0.5 )
458          bxr = bxl + blx
459
460          IF ( building_wall_south == 9999999.9 )  THEN
461             building_wall_south = ( ny + 1 - bly ) / 2 * dy
462          ENDIF
463          bys = INT ( building_wall_south / dy + 0.5 )
464          byn = bys + bly
465
466          WRITE ( io, 271 )  building_length_x, building_length_y, &
467                             building_height, bxl, bxr, bys, byn
468
469       CASE ( 'single_street_canyon' )
470          ch  = NINT( canyon_height / dz )
471          IF ( canyon_width_x /= 9999999.9 )  THEN
472!
473!--          Street canyon in y direction
474             cwx = NINT( canyon_width_x / dx )
475             IF ( canyon_wall_left == 9999999.9 )  THEN
476                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
477             ENDIF
478             cxl = NINT( canyon_wall_left / dx )
479             cxr = cxl + cwx
480             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
481
482          ELSEIF ( canyon_width_y /= 9999999.9 )  THEN
483!
484!--          Street canyon in x direction
485             cwy = NINT( canyon_width_y / dy )
486             IF ( canyon_wall_south == 9999999.9 )  THEN
487                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
488             ENDIF
489             cys = NINT( canyon_wall_south / dy )
490             cyn = cys + cwy
491             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
492          ENDIF
493
494    END SELECT
495
496    IF ( TRIM( topography ) /= 'flat' )  THEN
497       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
498          IF ( TRIM( topography ) == 'single_building' .OR.  &
499               TRIM( topography ) == 'single_street_canyon' )  THEN
500             WRITE ( io, 278 )
501          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
502             WRITE ( io, 279 )
503          ENDIF
504       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
505          WRITE ( io, 278 )
506       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
507          WRITE ( io, 279 )
508       ENDIF
509    ENDIF
510
511    IF ( plant_canopy ) THEN
512
513       WRITE ( io, 280 ) canopy_mode, pch_index, drag_coefficient
514       IF ( passive_scalar ) THEN
515          WRITE ( io, 281 ) scalar_exchange_coefficient,   &
516                            leaf_surface_concentration
517       ENDIF
518
519!
520!--    Heat flux at the top of vegetation
521       WRITE ( io, 282 ) cthf
522
523!
524!--    Leaf area density profile
525!--    Building output strings, starting with surface value
526       WRITE ( learde, '(F6.4)' )  lad_surface
527       gradients = '------'
528       slices = '     0'
529       coordinates = '   0.0'
530       i = 1
531       DO  WHILE ( lad_vertical_gradient_level_ind(i) /= -9999 )
532
533          WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
534          learde = TRIM( learde ) // ' ' // TRIM( coor_chr )
535
536          WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
537          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
538
539          WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
540          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
541
542          WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
543          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
544
545          i = i + 1
546       ENDDO
547
548       WRITE ( io, 283 )  TRIM( coordinates ), TRIM( learde ), &
549                          TRIM( gradients ), TRIM( slices )
550
551    ENDIF
552
553!
554!-- Boundary conditions
555    IF ( ibc_p_b == 0 )  THEN
556       runten = 'p(0)     = 0      |'
557    ELSEIF ( ibc_p_b == 1 )  THEN
558       runten = 'p(0)     = p(1)   |'
559    ELSE
560       runten = 'p(0)     = p(1) +R|'
561    ENDIF
562    IF ( ibc_p_t == 0 )  THEN
563       roben  = 'p(nzt+1) = 0      |'
564    ELSE
565       roben  = 'p(nzt+1) = p(nzt) |'
566    ENDIF
567
568    IF ( ibc_uv_b == 0 )  THEN
569       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
570    ELSE
571       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
572    ENDIF
573    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
574       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
575    ELSEIF ( ibc_uv_t == 0 )  THEN
576       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
577    ELSE
578       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
579    ENDIF
580
581    IF ( ibc_pt_b == 0 )  THEN
582       runten = TRIM( runten ) // ' pt(0)   = pt_surface'
583    ELSEIF ( ibc_pt_b == 1 )  THEN
584       runten = TRIM( runten ) // ' pt(0)   = pt(1)'
585    ELSEIF ( ibc_pt_b == 2 )  THEN
586       runten = TRIM( runten ) // ' pt(0) = from coupled model'
587    ENDIF
588    IF ( ibc_pt_t == 0 )  THEN
589       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
590    ELSEIF( ibc_pt_t == 1 )  THEN
591       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
592    ELSEIF( ibc_pt_t == 2 )  THEN
593       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
594
595    ENDIF
596
597    WRITE ( io, 300 )  runten, roben
598
599    IF ( .NOT. constant_diffusion )  THEN
600       IF ( ibc_e_b == 1 )  THEN
601          runten = 'e(0)     = e(1)'
602       ELSE
603          runten = 'e(0)     = e(1) = (u*/0.1)**2'
604       ENDIF
605       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
606
607       WRITE ( io, 301 )  'e', runten, roben       
608
609    ENDIF
610
611    IF ( ocean )  THEN
612       runten = 'sa(0)    = sa(1)'
613       IF ( ibc_sa_t == 0 )  THEN
614          roben =  'sa(nzt+1) = sa_surface'
615       ELSE
616          roben =  'sa(nzt+1) = sa(nzt)'
617       ENDIF
618       WRITE ( io, 301 ) 'sa', runten, roben
619    ENDIF
620
621    IF ( humidity )  THEN
622       IF ( ibc_q_b == 0 )  THEN
623          runten = 'q(0)     = q_surface'
624       ELSE
625          runten = 'q(0)     = q(1)'
626       ENDIF
627       IF ( ibc_q_t == 0 )  THEN
628          roben =  'q(nzt)   = q_top'
629       ELSE
630          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
631       ENDIF
632       WRITE ( io, 301 ) 'q', runten, roben
633    ENDIF
634
635    IF ( passive_scalar )  THEN
636       IF ( ibc_q_b == 0 )  THEN
637          runten = 's(0)     = s_surface'
638       ELSE
639          runten = 's(0)     = s(1)'
640       ENDIF
641       IF ( ibc_q_t == 0 )  THEN
642          roben =  's(nzt)   = s_top'
643       ELSE
644          roben =  's(nzt)   = s(nzt-1) + ds/dz'
645       ENDIF
646       WRITE ( io, 301 ) 's', runten, roben
647    ENDIF
648
649    IF ( use_surface_fluxes )  THEN
650       WRITE ( io, 303 )
651       IF ( constant_heatflux )  THEN
652          WRITE ( io, 306 )  surface_heatflux
653          IF ( random_heatflux )  WRITE ( io, 307 )
654       ENDIF
655       IF ( humidity  .AND.  constant_waterflux )  THEN
656          WRITE ( io, 311 ) surface_waterflux
657       ENDIF
658       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
659          WRITE ( io, 313 ) surface_waterflux
660       ENDIF
661    ENDIF
662
663    IF ( use_top_fluxes )  THEN
664       WRITE ( io, 304 )
665       IF ( coupling_mode == 'uncoupled' )  THEN
666          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
667          IF ( constant_top_heatflux )  THEN
668             WRITE ( io, 306 )  top_heatflux
669          ENDIF
670       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
671          WRITE ( io, 316 )
672       ENDIF
673       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
674          WRITE ( io, 309 )  top_salinityflux
675       ENDIF
676       IF ( humidity  .OR.  passive_scalar )  THEN
677          WRITE ( io, 315 )
678       ENDIF
679    ENDIF
680
681    IF ( prandtl_layer )  THEN
682       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length, &
683                          z0h_factor*roughness_length, kappa, &
684                          rif_min, rif_max
685       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
686       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
687          WRITE ( io, 312 )
688       ENDIF
689       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
690          WRITE ( io, 314 )
691       ENDIF
692    ELSE
693       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
694          WRITE ( io, 310 )  rif_min, rif_max
695       ENDIF
696    ENDIF
697
698    WRITE ( io, 317 )  bc_lr, bc_ns
699    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
700       WRITE ( io, 318 )  pt_damping_width, pt_damping_factor       
701       IF ( turbulent_inflow )  THEN
702          WRITE ( io, 319 )  recycling_width, recycling_plane, &
703                             inflow_damping_height, inflow_damping_width
704       ENDIF
705    ENDIF
706
707!
708!-- Listing of 1D-profiles
709    WRITE ( io, 325 )  dt_dopr_listing
710    IF ( averaging_interval_pr /= 0.0 )  THEN
711       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
712    ENDIF
713
714!
715!-- DATA output
716    WRITE ( io, 330 )
717    IF ( averaging_interval_pr /= 0.0 )  THEN
718       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
719    ENDIF
720
721!
722!-- 1D-profiles
723    dopr_chr = 'Profile:'
724    IF ( dopr_n /= 0 )  THEN
725       WRITE ( io, 331 )
726
727       output_format = ''
728       IF ( netcdf_output )  THEN
729          IF ( netcdf_data_format == 1 )  THEN
730             output_format = 'NetCDF classic'
731          ELSE
732             output_format = 'NetCDF 64bit offset'
733          ENDIF
734       ENDIF
735       WRITE ( io, 344 )  output_format
736
737       DO  i = 1, dopr_n
738          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
739          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
740             WRITE ( io, 332 )  dopr_chr
741             dopr_chr = '       :'
742          ENDIF
743       ENDDO
744
745       IF ( dopr_chr /= '' )  THEN
746          WRITE ( io, 332 )  dopr_chr
747       ENDIF
748       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
749       IF ( skip_time_dopr /= 0.0 )  WRITE ( io, 339 )  skip_time_dopr
750    ENDIF
751
752!
753!-- 2D-arrays
754    DO  av = 0, 1
755
756       i = 1
757       do2d_xy = ''
758       do2d_xz = ''
759       do2d_yz = ''
760       DO  WHILE ( do2d(av,i) /= ' ' )
761
762          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
763          do2d_mode = do2d(av,i)(l-1:l)
764
765          SELECT CASE ( do2d_mode )
766             CASE ( 'xy' )
767                ll = LEN_TRIM( do2d_xy )
768                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
769             CASE ( 'xz' )
770                ll = LEN_TRIM( do2d_xz )
771                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
772             CASE ( 'yz' )
773                ll = LEN_TRIM( do2d_yz )
774                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
775          END SELECT
776
777          i = i + 1
778
779       ENDDO
780
781       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
782              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
783              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) )  .AND. &
784            ( netcdf_output  .OR.  iso2d_output ) )  THEN
785
786          IF (  av == 0 )  THEN
787             WRITE ( io, 334 )  ''
788          ELSE
789             WRITE ( io, 334 )  '(time-averaged)'
790          ENDIF
791
792          IF ( do2d_at_begin )  THEN
793             begin_chr = 'and at the start'
794          ELSE
795             begin_chr = ''
796          ENDIF
797
798          output_format = ''
799          IF ( netcdf_output )  THEN
800             IF ( netcdf_data_format == 1 )  THEN
801                output_format = 'NetCDF classic'
802             ELSEIF ( netcdf_data_format == 2 )  THEN
803                output_format = 'NetCDF 64bit offset'
804             ELSEIF ( netcdf_data_format == 3 )  THEN
805                output_format = 'NetCDF4/HDF5'
806             ELSEIF ( netcdf_data_format == 4 )  THEN
807                output_format = 'NetCDF4/HDF5 clasic'
808             ENDIF
809          ENDIF
810          IF ( iso2d_output )  THEN
811             IF ( netcdf_output )  THEN
812                output_format = TRIM( output_format ) // ' and iso2d'
813             ELSE
814                output_format = 'iso2d'
815             ENDIF
816          ENDIF
817          WRITE ( io, 344 )  output_format
818
819          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
820             i = 1
821             slices = '/'
822             coordinates = '/'
823!
824!--          Building strings with index and coordinate informations of the
825!--          slices
826             DO  WHILE ( section(i,1) /= -9999 )
827
828                WRITE (section_chr,'(I5)')  section(i,1)
829                section_chr = ADJUSTL( section_chr )
830                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
831
832                IF ( section(i,1) == -1 )  THEN
833                   WRITE (coor_chr,'(F10.1)')  -1.0
834                ELSE
835                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
836                ENDIF
837                coor_chr = ADJUSTL( coor_chr )
838                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
839
840                i = i + 1
841             ENDDO
842             IF ( av == 0 )  THEN
843                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
844                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
845                                   TRIM( coordinates )
846                IF ( skip_time_do2d_xy /= 0.0 )  THEN
847                   WRITE ( io, 339 )  skip_time_do2d_xy
848                ENDIF
849             ELSE
850                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
851                                   TRIM( begin_chr ), averaging_interval, &
852                                   dt_averaging_input, 'k', TRIM( slices ), &
853                                   TRIM( coordinates )
854                IF ( skip_time_data_output_av /= 0.0 )  THEN
855                   WRITE ( io, 339 )  skip_time_data_output_av
856                ENDIF
857             ENDIF
858
859          ENDIF
860
861          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
862             i = 1
863             slices = '/'
864             coordinates = '/'
865!
866!--          Building strings with index and coordinate informations of the
867!--          slices
868             DO  WHILE ( section(i,2) /= -9999 )
869
870                WRITE (section_chr,'(I5)')  section(i,2)
871                section_chr = ADJUSTL( section_chr )
872                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
873
874                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
875                coor_chr = ADJUSTL( coor_chr )
876                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
877
878                i = i + 1
879             ENDDO
880             IF ( av == 0 )  THEN
881                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
882                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
883                                   TRIM( coordinates )
884                IF ( skip_time_do2d_xz /= 0.0 )  THEN
885                   WRITE ( io, 339 )  skip_time_do2d_xz
886                ENDIF
887             ELSE
888                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
889                                   TRIM( begin_chr ), averaging_interval, &
890                                   dt_averaging_input, 'j', TRIM( slices ), &
891                                   TRIM( coordinates )
892                IF ( skip_time_data_output_av /= 0.0 )  THEN
893                   WRITE ( io, 339 )  skip_time_data_output_av
894                ENDIF
895             ENDIF
896          ENDIF
897
898          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
899             i = 1
900             slices = '/'
901             coordinates = '/'
902!
903!--          Building strings with index and coordinate informations of the
904!--          slices
905             DO  WHILE ( section(i,3) /= -9999 )
906
907                WRITE (section_chr,'(I5)')  section(i,3)
908                section_chr = ADJUSTL( section_chr )
909                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
910
911                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
912                coor_chr = ADJUSTL( coor_chr )
913                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
914
915                i = i + 1
916             ENDDO
917             IF ( av == 0 )  THEN
918                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
919                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
920                                   TRIM( coordinates )
921                IF ( skip_time_do2d_yz /= 0.0 )  THEN
922                   WRITE ( io, 339 )  skip_time_do2d_yz
923                ENDIF
924             ELSE
925                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
926                                   TRIM( begin_chr ), averaging_interval, &
927                                   dt_averaging_input, 'i', TRIM( slices ), &
928                                   TRIM( coordinates )
929                IF ( skip_time_data_output_av /= 0.0 )  THEN
930                   WRITE ( io, 339 )  skip_time_data_output_av
931                ENDIF
932             ENDIF
933          ENDIF
934
935       ENDIF
936
937    ENDDO
938
939!
940!-- 3d-arrays
941    DO  av = 0, 1
942
943       i = 1
944       do3d_chr = ''
945       DO  WHILE ( do3d(av,i) /= ' ' )
946
947          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
948          i = i + 1
949
950       ENDDO
951
952       IF ( do3d_chr /= '' )  THEN
953          IF ( av == 0 )  THEN
954             WRITE ( io, 336 )  ''
955          ELSE
956             WRITE ( io, 336 )  '(time-averaged)'
957          ENDIF
958
959          output_format = ''
960          IF ( netcdf_output )  THEN
961             IF ( netcdf_data_format == 1 )  THEN
962                output_format = 'NetCDF classic'
963             ELSEIF ( netcdf_data_format == 2 )  THEN
964                output_format = 'NetCDF 64bit offset'
965             ELSEIF ( netcdf_data_format == 3 )  THEN
966                output_format = 'NetCDF4/HDF5'
967             ELSEIF ( netcdf_data_format == 4 )  THEN
968                output_format = 'NetCDF4/HDF5 clasic'
969             ENDIF
970          ENDIF
971          IF ( avs_output )  THEN
972             IF ( netcdf_output )  THEN
973                output_format = TRIM( output_format ) // ' and avs'
974             ELSE
975                output_format = 'avs'
976             ENDIF
977          ENDIF
978          WRITE ( io, 344 )  output_format
979
980          IF ( do3d_at_begin )  THEN
981             begin_chr = 'and at the start'
982          ELSE
983             begin_chr = ''
984          ENDIF
985          IF ( av == 0 )  THEN
986             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
987                                zu(nz_do3d), nz_do3d
988          ELSE
989             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
990                                TRIM( begin_chr ), averaging_interval, &
991                                dt_averaging_input, zu(nz_do3d), nz_do3d
992          ENDIF
993
994          IF ( do3d_compress )  THEN
995             do3d_chr = ''
996             i = 1
997             DO WHILE ( do3d(av,i) /= ' ' )
998
999                SELECT CASE ( do3d(av,i) )
1000                   CASE ( 'u' )
1001                      j = 1
1002                   CASE ( 'v' )
1003                      j = 2
1004                   CASE ( 'w' )
1005                      j = 3
1006                   CASE ( 'p' )
1007                      j = 4
1008                   CASE ( 'pt' )
1009                      j = 5
1010                END SELECT
1011                WRITE ( prec, '(I1)' )  plot_3d_precision(j)%precision
1012                do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // &
1013                           ':' // prec // ','
1014                i = i + 1
1015
1016             ENDDO
1017             WRITE ( io, 338 )  do3d_chr
1018
1019          ENDIF
1020
1021          IF ( av == 0 )  THEN
1022             IF ( skip_time_do3d /= 0.0 )  THEN
1023                WRITE ( io, 339 )  skip_time_do3d
1024             ENDIF
1025          ELSE
1026             IF ( skip_time_data_output_av /= 0.0 )  THEN
1027                WRITE ( io, 339 )  skip_time_data_output_av
1028             ENDIF
1029          ENDIF
1030
1031       ENDIF
1032
1033    ENDDO
1034
1035!
1036!-- masked arrays
1037    IF ( masks > 0 )  WRITE ( io, 345 )  &
1038         mask_scale_x, mask_scale_y, mask_scale_z
1039    DO  mid = 1, masks
1040       DO  av = 0, 1
1041
1042          i = 1
1043          domask_chr = ''
1044          DO  WHILE ( domask(mid,av,i) /= ' ' )
1045             domask_chr = TRIM( domask_chr ) // ' ' //  &
1046                          TRIM( domask(mid,av,i) ) // ','
1047             i = i + 1
1048          ENDDO
1049
1050          IF ( domask_chr /= '' )  THEN
1051             IF ( av == 0 )  THEN
1052                WRITE ( io, 346 )  '', mid
1053             ELSE
1054                WRITE ( io, 346 )  ' (time-averaged)', mid
1055             ENDIF
1056
1057             output_format = ''
1058             IF ( netcdf_output )  THEN
1059                IF ( netcdf_data_format == 1 )  THEN
1060                   output_format = 'NetCDF classic'
1061                ELSEIF ( netcdf_data_format == 2 )  THEN
1062                   output_format = 'NetCDF 64bit offset'
1063                ELSEIF ( netcdf_data_format == 3 )  THEN
1064                   output_format = 'NetCDF4/HDF5'
1065                ELSEIF ( netcdf_data_format == 4 )  THEN
1066                   output_format = 'NetCDF4/HDF5 clasic'
1067                ENDIF
1068             ENDIF
1069             WRITE ( io, 344 )  output_format
1070
1071             IF ( av == 0 )  THEN
1072                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1073             ELSE
1074                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1075                                   averaging_interval, dt_averaging_input
1076             ENDIF
1077
1078             IF ( av == 0 )  THEN
1079                IF ( skip_time_domask(mid) /= 0.0 )  THEN
1080                   WRITE ( io, 339 )  skip_time_domask(mid)
1081                ENDIF
1082             ELSE
1083                IF ( skip_time_data_output_av /= 0.0 )  THEN
1084                   WRITE ( io, 339 )  skip_time_data_output_av
1085                ENDIF
1086             ENDIF
1087!
1088!--          output locations
1089             DO  dim = 1, 3
1090                IF ( mask(mid,dim,1) >= 0.0 )  THEN
1091                   count = 0
1092                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0 )
1093                      count = count + 1
1094                   ENDDO
1095                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1096                                      mask(mid,dim,:count)
1097                ELSEIF ( mask_loop(mid,dim,1) < 0.0 .AND.  &
1098                         mask_loop(mid,dim,2) < 0.0 .AND.  &
1099                         mask_loop(mid,dim,3) == 0.0 )  THEN
1100                   WRITE ( io, 350 )  dir(dim), dir(dim)
1101                ELSEIF ( mask_loop(mid,dim,3) == 0.0 )  THEN
1102                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1103                                      mask_loop(mid,dim,1:2)
1104                ELSE
1105                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1106                                      mask_loop(mid,dim,1:3)
1107                ENDIF
1108             ENDDO
1109          ENDIF
1110
1111       ENDDO
1112    ENDDO
1113
1114!
1115!-- Timeseries
1116    IF ( dt_dots /= 9999999.9 )  THEN
1117       WRITE ( io, 340 )
1118
1119       output_format = ''
1120       IF ( netcdf_output )  THEN
1121          IF ( netcdf_data_format == 1 )  THEN
1122             output_format = 'NetCDF classic'
1123          ELSE
1124             output_format = 'NetCDF 64bit offset'
1125          ENDIF
1126       ENDIF
1127       WRITE ( io, 344 )  output_format
1128       WRITE ( io, 341 )  dt_dots
1129    ENDIF
1130
1131#if defined( __dvrp_graphics )
1132!
1133!-- Dvrp-output
1134    IF ( dt_dvrp /= 9999999.9 )  THEN
1135       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1136                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1137       i = 1
1138       l = 0
1139       m = 0
1140       DO WHILE ( mode_dvrp(i) /= ' ' )
1141          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
1142             READ ( mode_dvrp(i), '(10X,I2)' )  j
1143             l = l + 1
1144             IF ( do3d(0,j) /= ' ' )  THEN
1145                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1146                                   isosurface_color(:,l)
1147             ENDIF
1148          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
1149             READ ( mode_dvrp(i), '(6X,I2)' )  j
1150             m = m + 1
1151             IF ( do2d(0,j) /= ' ' )  THEN
1152                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1153                                   slicer_range_limits_dvrp(:,m)
1154             ENDIF
1155          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
1156             WRITE ( io, 363 )  dvrp_psize
1157             IF ( particle_dvrpsize /= 'none' )  THEN
1158                WRITE ( io, 364 )  'size', TRIM( particle_dvrpsize ), &
1159                                   dvrpsize_interval
1160             ENDIF
1161             IF ( particle_color /= 'none' )  THEN
1162                WRITE ( io, 364 )  'color', TRIM( particle_color ), &
1163                                   color_interval
1164             ENDIF
1165          ENDIF
1166          i = i + 1
1167       ENDDO
1168
1169       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1170                          superelevation_y, superelevation, clip_dvrp_l, &
1171                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1172
1173       IF ( TRIM( topography ) /= 'flat' )  THEN
1174          WRITE ( io, 366 )  topography_color
1175          IF ( cluster_size > 1 )  THEN
1176             WRITE ( io, 367 )  cluster_size
1177          ENDIF
1178       ENDIF
1179
1180    ENDIF
1181#endif
1182
1183#if defined( __spectra )
1184!
1185!-- Spectra output
1186    IF ( dt_dosp /= 9999999.9 ) THEN
1187       WRITE ( io, 370 )
1188
1189       output_format = ''
1190       IF ( netcdf_output )  THEN
1191          IF ( netcdf_data_format == 1 )  THEN
1192             output_format = 'NetCDF classic'
1193          ELSE
1194             output_format = 'NetCDF 64bit offset'
1195          ENDIF
1196       ENDIF
1197       WRITE ( io, 344 )  output_format
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 ),  &
1202                          ( comp_spectra_level(i), i = 1,100 ), &
1203                          ( plot_spectra_level(i), i = 1,100 ), &
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
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
1228
1229!
1230!-- Cloud physics parameters
1231    IF ( cloud_physics ) THEN
1232       WRITE ( io, 415 )
1233       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
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     
1245       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
1246       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
1247
1248       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
1249       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1250
1251       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
1252       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1253
1254       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
1255       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1256
1257       IF ( i == 10 )  THEN
1258          EXIT
1259       ELSE
1260          i = i + 1
1261       ENDIF
1262
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
1277       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
1278       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
1279
1280       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
1281       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1282
1283       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
1284       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1285
1286       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
1287       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1288
1289       IF ( i == 10 )  THEN
1290          EXIT
1291       ELSE
1292          i = i + 1
1293       ENDIF
1294 
1295    ENDDO
1296
1297    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
1298                       TRIM( gradients ), TRIM( slices )
1299
1300!
1301!-- Initial wind profiles
1302    IF ( u_profile(1) /= 9999999.9 )  WRITE ( io, 427 )
1303
1304!
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
1314       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1315       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1316
1317       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1318       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1319
1320       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1321       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1322
1323       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1324       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1325
1326       IF ( i == 10 )  THEN
1327          EXIT
1328       ELSE
1329          i = i + 1
1330       ENDIF
1331
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
1340    IF ( humidity  .OR.  passive_scalar )  THEN
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
1360          IF ( i == 10 )  THEN
1361             EXIT
1362          ELSE
1363             i = i + 1
1364          ENDIF
1365
1366       ENDDO
1367
1368       IF ( humidity )  THEN
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!
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
1400          IF ( i == 10 )  THEN
1401             EXIT
1402          ELSE
1403             i = i + 1
1404          ENDIF
1405
1406       ENDDO
1407
1408       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1409                          TRIM( gradients ), TRIM( slices )
1410    ENDIF
1411
1412!
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
1421       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
1422
1423          WRITE (coor_chr,'(E10.2,7X)')  &
1424                                w_subs(subs_vertical_gradient_level_i(i))
1425          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1426
1427          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
1428          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1429
1430          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
1431          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1432
1433          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
1434          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1435
1436          IF ( i == 10 )  THEN
1437             EXIT
1438          ELSE
1439             i = i + 1
1440          ENDIF
1441
1442       ENDDO
1443
1444       WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
1445                          TRIM( gradients ), TRIM( slices )
1446    ENDIF
1447
1448!
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 )
1460       IF ( collision_kernel /= 'none' )  THEN
1461          WRITE ( io, 435 )  TRIM( collision_kernel )
1462          IF ( collision_kernel(6:9) == 'fast' )  THEN
1463             WRITE ( io, 436 )  radius_classes, dissipation_classes
1464          ENDIF
1465       ELSE
1466          WRITE ( io, 437 )
1467       ENDIF
1468    ENDIF
1469
1470!
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
1482       IF ( e_init > 0.0 )  WRITE ( io, 455 )  e_init
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
1494       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
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
1504    IF ( humidity  .AND.  q_surface_initial_change /= 0.0 )  THEN
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
1511    IF ( particle_advection )  THEN
1512!
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, &
1516                          end_time_prel, dt_sort_particles
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
1521       IF ( use_particle_tails  .AND.  maximum_number_of_tailpoints /= 0 )  THEN
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
1528       ENDIF
1529       IF ( dt_write_particle_data /= 9999999.9 )  THEN
1530          WRITE ( io, 485 )  dt_write_particle_data
1531          output_format = ''
1532          IF ( netcdf_output )  THEN
1533             IF ( netcdf_data_format > 1 )  THEN
1534                output_format = 'netcdf (64 bit offset) and binary'
1535             ELSE
1536                output_format = 'netcdf and binary'
1537             ENDIF
1538          ELSE
1539             output_format = 'binary'
1540          ENDIF
1541          WRITE ( io, 344 )  output_format
1542       ENDIF
1543       IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1544       IF ( write_particle_statistics )  WRITE ( io, 486 )
1545
1546       WRITE ( io, 487 )  number_of_particle_groups
1547
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 )
1552          ELSE
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
1559          ENDIF
1560          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1561                             pdx(i), pdy(i), pdz(i)
1562          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
1563       ENDDO
1564
1565    ENDIF
1566
1567
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
1586    CALL local_flush( io )
1587
1588!
1589!-- Here the FORMATs start
1590
1591 99 FORMAT (1X,78('-'))
1592100 FORMAT (/1X,'***************************',9X,42('-')/        &
1593            1X,'* ',A,' *',9X,A/                               &
1594            1X,'***************************',9X,42('-'))
1595101 FORMAT (37X,'coupled run using MPI-',I1,': ',A/ &
1596            37X,42('-'))
1597102 FORMAT (/' Date:              ',A8,9X,'Run:       ',A20/      &
1598            ' Time:              ',A8,9X,'Run-No.:   ',I2.2/     &
1599            ' Run on host:     ',A10)
1600#if defined( __parallel )
1601103 FORMAT (' Number of PEs:',8X,I5,9X,'Processor grid (x,y): (',I3,',',I3, &
1602              ')',1X,A)
1603104 FORMAT (' Number of PEs:',8X,I5,9X,'Tasks:',I4,'   threads per task:',I4/ &
1604              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
1605105 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1606106 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
1607            37X,'because the job is running on an SMP-cluster')
1608107 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
1609108 FORMAT (37X,'Max. # of parallel I/O streams is ',I5)
1610#endif
1611110 FORMAT (/' Numerical Schemes:'/ &
1612             ' -----------------'/)
1613111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1614112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1615            '     Iterations (initial/other): ',I3,'/',I3,'  omega = ',F5.3)
1616113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1617                  ' or Upstream')
1618116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1619                  ' or Upstream')
1620118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1621119 FORMAT (' --> Galilei-Transform applied to horizontal advection', &
1622            '     Translation velocity = ',A/ &
1623            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1624122 FORMAT (' --> Time differencing scheme: ',A)
1625123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
1626            '     maximum damping coefficient: ',F5.3, ' 1/s')
1627129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1628130 FORMAT (' --> Additional prognostic equation for the total water content')
1629131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
1630                  F6.2, ' K assumed')
1631132 FORMAT ('     Parameterization of long-wave radiation processes via'/ &
1632            '     effective emissivity scheme')
1633133 FORMAT ('     Precipitation parameterization via Kessler-Scheme')
1634134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1635135 FORMAT (' --> Solve perturbation pressure via multigrid method (', &
1636                  A,'-cycle)'/ &
1637            '     number of grid levels:                   ',I2/ &
1638            '     Gauss-Seidel red/black iterations:       ',I2)
1639136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1640                  I3,')')
1641137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1642            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1643                  I3,')'/ &
1644            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1645                  I3,')')
1646138 FORMAT ('     Using hybrid version for 1d-domain-decomposition')
1647139 FORMAT (' --> Loop optimization method: ',A)
1648140 FORMAT ('     maximum residual allowed:                ',E10.3)
1649141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1650142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1651                  'step')
1652143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1653                  'kinetic energy')
1654144 FORMAT ('     masking method is used')
1655150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
1656                  'conserved'/ &
1657            '     using the ',A,' mode')
1658151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
1659152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
1660           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
1661           /'     starting from dp_level_b =', F8.3, 'm', A /)
1662153 FORMAT (' --> Large-scale vertical motion is used in the ', &
1663                  'prognostic equation for')
1664154 FORMAT ('     the potential temperature')
1665200 FORMAT (//' Run time and time step information:'/ &
1666             ' ----------------------------------'/)
1667201 FORMAT ( ' Timestep:          variable     maximum value: ',F6.3,' s', &
1668             '    CFL-factor: ',F4.2)
1669202 FORMAT ( ' Timestep:       dt = ',F6.3,' s'/)
1670203 FORMAT ( ' Start time:       ',F9.3,' s'/ &
1671             ' End time:         ',F9.3,' s')
1672204 FORMAT ( A,F9.3,' s')
1673205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
1674206 FORMAT (/' Time reached:     ',F9.3,' s'/ &
1675             ' CPU-time used:    ',F9.3,' s     per timestep:               ', &
1676               '  ',F9.3,' s'/                                                 &
1677             '                                   per second of simulated tim', &
1678               'e: ',F9.3,' s')
1679207 FORMAT ( A/' Coupling start time:',F9.3,' s')
1680250 FORMAT (//' Computational grid and domain size:'/ &
1681              ' ----------------------------------'// &
1682              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1683              ' m    dz =    ',F7.3,' m'/ &
1684              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1685              ' m  z(u) = ',F10.3,' m'/)
1686252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
1687              ' factor: ',F5.3/ &
1688            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1689254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1690            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1691260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1692             ' degrees')
1693270 FORMAT (//' Topography informations:'/ &
1694              ' -----------------------'// &
1695              1X,'Topography: ',A)
1696271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1697              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1698                ' / ',I4)
1699272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
1700              ' direction' / &
1701              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
1702              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
1703278 FORMAT (' Topography grid definition convention:'/ &
1704            ' cell edge (staggered grid points'/  &
1705            ' (u in x-direction, v in y-direction))' /)
1706279 FORMAT (' Topography grid definition convention:'/ &
1707            ' cell center (scalar grid points)' /)
1708280 FORMAT (//' Vegetation canopy (drag) model:'/ &
1709              ' ------------------------------'// &
1710              ' Canopy mode: ', A / &
1711              ' Canopy top: ',I4 / &
1712              ' Leaf drag coefficient: ',F6.2 /)
1713281 FORMAT (/ ' Scalar_exchange_coefficient: ',F6.2 / &
1714              ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /)
1715282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s')
1716283 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
1717              ' Height:              ',A,'  m'/ &
1718              ' Leaf area density:   ',A,'  m**2/m**3'/ &
1719              ' Gradient:            ',A,'  m**2/m**4'/ &
1720              ' Gridpoint:           ',A)
1721               
1722300 FORMAT (//' Boundary conditions:'/ &
1723             ' -------------------'// &
1724             '                     p                    uv             ', &
1725             '                   pt'// &
1726             ' B. bound.: ',A/ &
1727             ' T. bound.: ',A)
1728301 FORMAT (/'                     ',A// &
1729             ' B. bound.: ',A/ &
1730             ' T. bound.: ',A)
1731303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1732304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1733305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
1734               'computational u,v-level:'// &
1735             '       zp = ',F6.2,' m   z0 = ',F6.4,' m   z0h = ',F7.5,&
1736             ' m   kappa = ',F4.2/ &
1737             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1738306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
1739307 FORMAT ('       Heatflux has a random normal distribution')
1740308 FORMAT ('       Predefined surface temperature')
1741309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
1742310 FORMAT (//'    1D-Model:'// &
1743             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1744311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1745312 FORMAT ('       Predefined surface humidity')
1746313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1747314 FORMAT ('       Predefined scalar value at the surface')
1748315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
1749316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1750                    'atmosphere model')
1751317 FORMAT (//' Lateral boundaries:'/ &
1752            '       left/right:  ',A/    &
1753            '       north/south: ',A)
1754318 FORMAT (/'       pt damping layer width = ',F7.2,' m, pt ', &
1755                    'damping factor = ',F6.4)
1756319 FORMAT ('       turbulence recycling at inflow switched on'/ &
1757            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1758            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
1759320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
1760            '                                          v: ',F9.6,' m**2/s**2')
1761325 FORMAT (//' List output:'/ &
1762             ' -----------'//  &
1763            '    1D-Profiles:'/    &
1764            '       Output every             ',F8.2,' s')
1765326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
1766            '       Averaging input every    ',F8.2,' s')
1767330 FORMAT (//' Data output:'/ &
1768             ' -----------'/)
1769331 FORMAT (/'    1D-Profiles:')
1770332 FORMAT (/'       ',A)
1771333 FORMAT ('       Output every             ',F8.2,' s',/ &
1772            '       Time averaged over       ',F8.2,' s'/ &
1773            '       Averaging input every    ',F8.2,' s')
1774334 FORMAT (/'    2D-Arrays',A,':')
1775335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1776            '       Output every             ',F8.2,' s  ',A/ &
1777            '       Cross sections at ',A1,' = ',A/ &
1778            '       scalar-coordinates:   ',A,' m'/)
1779336 FORMAT (/'    3D-Arrays',A,':')
1780337 FORMAT (/'       Arrays: ',A/ &
1781            '       Output every             ',F8.2,' s  ',A/ &
1782            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1783338 FORMAT ('       Compressed data output'/ &
1784            '       Decimal precision: ',A/)
1785339 FORMAT ('       No output during initial ',F8.2,' s')
1786340 FORMAT (/'    Time series:')
1787341 FORMAT ('       Output every             ',F8.2,' s'/)
1788342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1789            '       Output every             ',F8.2,' s  ',A/ &
1790            '       Time averaged over       ',F8.2,' s'/ &
1791            '       Averaging input every    ',F8.2,' s'/ &
1792            '       Cross sections at ',A1,' = ',A/ &
1793            '       scalar-coordinates:   ',A,' m'/)
1794343 FORMAT (/'       Arrays: ',A/ &
1795            '       Output every             ',F8.2,' s  ',A/ &
1796            '       Time averaged over       ',F8.2,' s'/ &
1797            '       Averaging input every    ',F8.2,' s'/ &
1798            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1799344 FORMAT ('       Output format: ',A/)
1800345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
1801            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
1802            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
1803            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
1804346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
1805347 FORMAT ('       Variables: ',A/ &
1806            '       Output every             ',F8.2,' s')
1807348 FORMAT ('       Variables: ',A/ &
1808            '       Output every             ',F8.2,' s'/ &
1809            '       Time averaged over       ',F8.2,' s'/ &
1810            '       Averaging input every    ',F8.2,' s')
1811349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1812            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
1813            13('       ',8(F8.2,',')/) )
1814350 FORMAT (/'       Output locations in ',A,'-direction: ', &
1815            'all gridpoints along ',A,'-direction (default).' )
1816351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1817            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
1818            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
1819#if defined( __dvrp_graphics )
1820360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
1821            '       Output every      ',F7.1,' s'/ &
1822            '       Output mode:      ',A/ &
1823            '       Host / User:      ',A,' / ',A/ &
1824            '       Directory:        ',A// &
1825            '       The sequence contains:')
1826361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
1827            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
1828362 FORMAT (/'       Slicer plane ',A/ &
1829            '       Slicer limits: [',F6.2,',',F6.2,']')
1830363 FORMAT (/'       Particles'/ &
1831            '          particle size:  ',F7.2,' m')
1832364 FORMAT ('          particle ',A,' controlled by "',A,'" with interval [', &
1833                       F6.2,',',F6.2,']')
1834365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
1835            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
1836                     ')'/ &
1837            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
1838            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
1839366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
1840367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
1841#endif
1842#if defined( __spectra )
1843370 FORMAT ('    Spectra:')
1844371 FORMAT ('       Output every ',F7.1,' s'/)
1845372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
1846            '       Directions: ', 10(A5,',')/                         &
1847            '       height levels  k = ', 20(I3,',')/                  &
1848            '                          ', 20(I3,',')/                  &
1849            '                          ', 20(I3,',')/                  &
1850            '                          ', 20(I3,',')/                  &
1851            '                          ', 19(I3,','),I3,'.'/           &
1852            '       height levels selected for standard plot:'/        &
1853            '                      k = ', 20(I3,',')/                  &
1854            '                          ', 20(I3,',')/                  &
1855            '                          ', 20(I3,',')/                  &
1856            '                          ', 20(I3,',')/                  &
1857            '                          ', 19(I3,','),I3,'.'/           &
1858            '       Time averaged over ', F7.1, ' s,' /                &
1859            '       Profiles for the time averaging are taken every ', &
1860                    F6.1,' s')
1861#endif
1862400 FORMAT (//' Physical quantities:'/ &
1863              ' -------------------'/)
1864410 FORMAT ('    Angular velocity    :   omega = ',E9.3,' rad/s'/  &
1865            '    Geograph. latitude  :   phi   = ',F4.1,' degr'/   &
1866            '    Coriolis parameter  :   f     = ',F9.6,' 1/s'/    &
1867            '                            f*    = ',F9.6,' 1/s')
1868411 FORMAT (/'    Gravity             :   g     = ',F4.1,' m/s**2')
1869412 FORMAT (/'    Reference density in buoyancy terms: ',F8.3,' kg/m**3')
1870413 FORMAT (/'    Reference temperature in buoyancy terms: ',F8.4,' K')
1871415 FORMAT (/'    Cloud physics parameters:'/ &
1872             '    ------------------------'/)
1873416 FORMAT ('        Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
1874            '        Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
1875            '        Density of air     :   rho_0 = ',F5.3,' kg/m**3'/  &
1876            '        Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
1877            '        Vapourization heat :   L_v   = ',E8.2,' J/kg')
1878420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1879            '       Height:        ',A,'  m'/ &
1880            '       Temperature:   ',A,'  K'/ &
1881            '       Gradient:      ',A,'  K/100m'/ &
1882            '       Gridpoint:     ',A)
1883421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1884            '       Height:      ',A,'  m'/ &
1885            '       Humidity:    ',A,'  kg/kg'/ &
1886            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1887            '       Gridpoint:   ',A)
1888422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1889            '       Height:                  ',A,'  m'/ &
1890            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1891            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1892            '       Gridpoint:               ',A)
1893423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1894            '       Height:      ',A,'  m'/ &
1895            '       ug:          ',A,'  m/s'/ &
1896            '       Gradient:    ',A,'  1/100s'/ &
1897            '       Gridpoint:   ',A)
1898424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1899            '       Height:      ',A,'  m'/ &
1900            '       vg:          ',A,'  m/s'/ &
1901            '       Gradient:    ',A,'  1/100s'/ &
1902            '       Gridpoint:   ',A)
1903425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
1904            '       Height:     ',A,'  m'/ &
1905            '       Salinity:   ',A,'  psu'/ &
1906            '       Gradient:   ',A,'  psu/100m'/ &
1907            '       Gridpoint:  ',A)
1908426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
1909            '       Height:      ',A,'  m'/ &
1910            '       w_subs:      ',A,'  m/s'/ &
1911            '       Gradient:    ',A,'  (m/s)/100m'/ &
1912            '       Gridpoint:   ',A)
1913427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
1914                  ' profiles')
1915430 FORMAT (//' Cloud physics quantities / methods:'/ &
1916              ' ----------------------------------'/)
1917431 FORMAT ('    Humidity is treated as purely passive scalar (no condensati', &
1918                 'on)')
1919432 FORMAT ('    Bulk scheme with liquid water potential temperature and'/ &
1920            '    total water content is used.'/ &
1921            '    Condensation is parameterized via 0% - or 100% scheme.')
1922433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
1923                 'icle model')
1924434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
1925                 ' droplets < 1.0E-6 m')
1926435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
1927436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
1928                    'are used'/ &
1929            '          number of radius classes:       ',I3,'    interval ', &
1930                       '[1.0E-6,2.0E-4] m'/ &
1931            '          number of dissipation classes:   ',I2,'    interval ', &
1932                       '[0,1000] cm**2/s**3')
1933437 FORMAT ('    Droplet collision is switched off')
1934450 FORMAT (//' LES / Turbulence quantities:'/ &
1935              ' ---------------------------'/)
1936451 FORMAT ('    Diffusion coefficients are constant:'/ &
1937            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1938453 FORMAT ('    Mixing length is limited to ',F4.2,' * z')
1939454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
1940455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
1941470 FORMAT (//' Actions during the simulation:'/ &
1942              ' -----------------------------'/)
1943471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1944            '    Disturbance amplitude           :     ',F4.2, ' m/s'/       &
1945            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1946            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
1947472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1948                 ' to i/j =',I4)
1949473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1950                 1X,F5.3, ' m**2/s**2')
1951474 FORMAT ('    Random number generator used    : ',A/)
1952475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1953                 'respectively, if'/ &
1954            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1955                 ' 3D-simulation'/)
1956476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1957                 'respectively, if the'/ &
1958            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1959                 ' the 3D-simulation'/)
1960477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1961                 'respectively, if the'/ &
1962            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1963                 ' the 3D-simulation'/)
1964480 FORMAT ('    Particles:'/ &
1965            '    ---------'// &
1966            '       Particle advection is active (switched on at t = ', F7.1, &
1967                    ' s)'/ &
1968            '       Start of new particle generations every  ',F6.1,' s'/ &
1969            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
1970            '                            bottom:     ', A, ' top:         ', A/&
1971            '       Maximum particle age:                 ',F9.1,' s'/ &
1972            '       Advection stopped at t = ',F9.1,' s'/ &
1973            '       Particles are sorted every ',F9.1,' s'/)
1974481 FORMAT ('       Particles have random start positions'/)
1975482 FORMAT ('          Particles are advected only horizontally'/)
1976483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
1977484 FORMAT ('            Number of tails of the total domain: ',I10/ &
1978            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
1979            '            Maximum age of the end of the tail:  ',F8.2,' s')
1980485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
1981486 FORMAT ('       Particle statistics are written on file'/)
1982487 FORMAT ('       Number of particle groups: ',I2/)
1983488 FORMAT ('       SGS velocity components are used for particle advection'/ &
1984            '          minimum timestep for advection: ', F7.5/)
1985489 FORMAT ('       Number of particles simultaneously released at each ', &
1986                    'point: ', I5/)
1987490 FORMAT ('       Particle group ',I2,':'/ &
1988            '          Particle radius: ',E10.3, 'm')
1989491 FORMAT ('          Particle inertia is activated'/ &
1990            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
1991492 FORMAT ('          Particles are advected only passively (no inertia)'/)
1992493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
1993            '                                         y:',F8.1,' - ',F8.1,' m'/&
1994            '                                         z:',F8.1,' - ',F8.1,' m'/&
1995            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
1996                       ' m  dz = ',F8.1,' m'/)
1997494 FORMAT ('       Output of particle time series in NetCDF format every ', &
1998                    F8.2,' s'/)
1999495 FORMAT ('       Number of particles in total domain: ',I10/)
2000500 FORMAT (//' 1D-Model parameters:'/                           &
2001              ' -------------------'//                           &
2002            '    Simulation time:                   ',F8.1,' s'/ &
2003            '    Run-controll output every:         ',F8.1,' s'/ &
2004            '    Vertical profile output every:     ',F8.1,' s'/ &
2005            '    Mixing length calculation:         ',A/         &
2006            '    Dissipation calculation:           ',A/)
2007502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
2008503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2009504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
2010
2011
2012 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.