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

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

revisions r1031 and r1034 documented

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