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

Last change on this file since 821 was 768, checked in by raasch, 13 years ago

last commit documented

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