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

Last change on this file since 448 was 448, checked in by raasch, 14 years ago

test print* removed in header

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