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

Last change on this file since 431 was 430, checked in by raasch, 15 years ago

bugfix concerning gradient_level indices in header

  • Property svn:keywords set to Id
File size: 75.9 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 430 2010-01-29 09:14:09Z 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       print*, 'i=',i,' level_u=',ug_vertical_gradient_level_ind(i)
1232       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
1233       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
1234
1235       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
1236       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1237
1238       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
1239       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1240
1241       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
1242       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1243
1244       IF ( i == 10 )  THEN
1245          EXIT
1246       ELSE
1247          i = i + 1
1248       ENDIF
1249
1250    ENDDO
1251
1252    WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
1253                       TRIM( gradients ), TRIM( slices )
1254
1255!-- Profile of the geostrophic wind (component vg)
1256!-- Building output strings
1257    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
1258    gradients = '------'
1259    slices = '     0'
1260    coordinates = '   0.0'
1261    i = 1
1262    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
1263
1264       print*, 'i=',i,' level_v=',vg_vertical_gradient_level_ind(i)
1265       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
1266       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
1267
1268       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
1269       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1270
1271       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
1272       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1273
1274       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
1275       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1276
1277       IF ( i == 10 )  THEN
1278          EXIT
1279       ELSE
1280          i = i + 1
1281       ENDIF
1282 
1283    ENDDO
1284
1285    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
1286                       TRIM( gradients ), TRIM( slices )
1287
1288!
1289!-- Initial temperature profile
1290!-- Building output strings, starting with surface temperature
1291    WRITE ( temperatures, '(F6.2)' )  pt_surface
1292    gradients = '------'
1293    slices = '     0'
1294    coordinates = '   0.0'
1295    i = 1
1296    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1297
1298       print*, 'i=',i,' level_pt=',pt_vertical_gradient_level_ind(i)
1299       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1300       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1301
1302       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1303       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1304
1305       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1306       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1307
1308       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1309       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1310
1311       IF ( i == 10 )  THEN
1312          EXIT
1313       ELSE
1314          i = i + 1
1315       ENDIF
1316
1317    ENDDO
1318
1319    WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1320                       TRIM( gradients ), TRIM( slices )
1321
1322!
1323!-- Initial humidity profile
1324!-- Building output strings, starting with surface humidity
1325    IF ( humidity  .OR.  passive_scalar )  THEN
1326       WRITE ( temperatures, '(E8.1)' )  q_surface
1327       gradients = '--------'
1328       slices = '       0'
1329       coordinates = '     0.0'
1330       i = 1
1331       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1332         
1333          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1334          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1335
1336          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1337          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1338         
1339          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1340          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1341         
1342          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1343          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1344
1345          IF ( i == 10 )  THEN
1346             EXIT
1347          ELSE
1348             i = i + 1
1349          ENDIF
1350
1351       ENDDO
1352
1353       IF ( humidity )  THEN
1354          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1355                             TRIM( gradients ), TRIM( slices )
1356       ELSE
1357          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1358                             TRIM( gradients ), TRIM( slices )
1359       ENDIF
1360    ENDIF
1361
1362!
1363!-- Initial salinity profile
1364!-- Building output strings, starting with surface salinity
1365    IF ( ocean )  THEN
1366       WRITE ( temperatures, '(F6.2)' )  sa_surface
1367       gradients = '------'
1368       slices = '     0'
1369       coordinates = '   0.0'
1370       i = 1
1371       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1372
1373          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1374          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1375
1376          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1377          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1378
1379          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1380          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1381
1382          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1383          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1384
1385          IF ( i == 10 )  THEN
1386             EXIT
1387          ELSE
1388             i = i + 1
1389          ENDIF
1390
1391       ENDDO
1392
1393       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1394                          TRIM( gradients ), TRIM( slices )
1395    ENDIF
1396
1397!
1398!-- Profile for the large scale vertial velocity
1399!-- Building output strings, starting with surface value
1400    IF ( large_scale_subsidence )  THEN
1401       temperatures = '   0.0'
1402       gradients = '------'
1403       slices = '     0'
1404       coordinates = '   0.0'
1405       i = 1
1406       DO  WHILE ( ws_vertical_gradient_level_ind(i) /= -9999 )
1407
1408          WRITE (coor_chr,'(E10.2,7X)')  &
1409                                w_subs(ws_vertical_gradient_level_ind(i))
1410          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1411
1412          WRITE (coor_chr,'(E10.2,7X)')  ws_vertical_gradient(i)
1413          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1414
1415          WRITE (coor_chr,'(I10,7X)')  ws_vertical_gradient_level_ind(i)
1416          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1417
1418          WRITE (coor_chr,'(F10.2,7X)')  ws_vertical_gradient_level(i)
1419          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1420
1421          IF ( i == 10 )  THEN
1422             EXIT
1423          ELSE
1424             i = i + 1
1425          ENDIF
1426
1427       ENDDO
1428
1429       WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
1430                          TRIM( gradients ), TRIM( slices )
1431    ENDIF
1432
1433!
1434!-- LES / turbulence parameters
1435    WRITE ( io, 450 )
1436
1437!--
1438! ... LES-constants used must still be added here
1439!--
1440    IF ( constant_diffusion )  THEN
1441       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1442                          prandtl_number
1443    ENDIF
1444    IF ( .NOT. constant_diffusion)  THEN
1445       IF ( e_init > 0.0 )  WRITE ( io, 455 )  e_init
1446       IF ( e_min > 0.0 )  WRITE ( io, 454 )  e_min
1447       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1448       IF ( adjust_mixing_length  .AND.  prandtl_layer )  WRITE ( io, 452 )
1449    ENDIF
1450
1451!
1452!-- Special actions during the run
1453    WRITE ( io, 470 )
1454    IF ( create_disturbances )  THEN
1455       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1456                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1457                          zu(disturbance_level_ind_t), disturbance_level_ind_t
1458       IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
1459          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1460       ELSE
1461          WRITE ( io, 473 )  disturbance_energy_limit
1462       ENDIF
1463       WRITE ( io, 474 )  TRIM( random_generator )
1464    ENDIF
1465    IF ( pt_surface_initial_change /= 0.0 )  THEN
1466       WRITE ( io, 475 )  pt_surface_initial_change
1467    ENDIF
1468    IF ( humidity  .AND.  q_surface_initial_change /= 0.0 )  THEN
1469       WRITE ( io, 476 )  q_surface_initial_change       
1470    ENDIF
1471    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0 )  THEN
1472       WRITE ( io, 477 )  q_surface_initial_change       
1473    ENDIF
1474
1475    IF ( particle_advection )  THEN
1476!
1477!--    Particle attributes
1478       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1479                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
1480                          end_time_prel, dt_sort_particles
1481       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1482       IF ( random_start_position )  WRITE ( io, 481 )
1483       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1484       WRITE ( io, 495 )  total_number_of_particles
1485       IF ( maximum_number_of_tailpoints /= 0 )  THEN
1486          WRITE ( io, 483 )  maximum_number_of_tailpoints
1487          IF ( minimum_tailpoint_distance /= 0 )  THEN
1488             WRITE ( io, 484 )  total_number_of_tails,      &
1489                                minimum_tailpoint_distance, &
1490                                maximum_tailpoint_age
1491          ENDIF
1492       ENDIF
1493       IF ( dt_write_particle_data /= 9999999.9 )  THEN
1494          WRITE ( io, 485 )  dt_write_particle_data
1495          output_format = ''
1496          IF ( netcdf_output )  THEN
1497             IF ( netcdf_64bit )  THEN
1498                output_format = 'netcdf (64 bit offset) and binary'
1499             ELSE
1500                output_format = 'netcdf and binary'
1501             ENDIF
1502          ELSE
1503             output_format = 'binary'
1504          ENDIF
1505          WRITE ( io, 344 )  output_format
1506       ENDIF
1507       IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1508       IF ( write_particle_statistics )  WRITE ( io, 486 )
1509
1510       WRITE ( io, 487 )  number_of_particle_groups
1511
1512       DO  i = 1, number_of_particle_groups
1513          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9 )  THEN
1514             WRITE ( io, 490 )  i, 0.0
1515             WRITE ( io, 492 )
1516          ELSE
1517             WRITE ( io, 490 )  i, radius(i)
1518             IF ( density_ratio(i) /= 0.0 )  THEN
1519                WRITE ( io, 491 )  density_ratio(i)
1520             ELSE
1521                WRITE ( io, 492 )
1522             ENDIF
1523          ENDIF
1524          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1525                             pdx(i), pdy(i), pdz(i)
1526          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
1527       ENDDO
1528
1529    ENDIF
1530
1531
1532!
1533!-- Parameters of 1D-model
1534    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1535       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1536                          mixing_length_1d, dissipation_1d
1537       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1538          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1539       ENDIF
1540    ENDIF
1541
1542!
1543!-- User-defined informations
1544    CALL user_header( io )
1545
1546    WRITE ( io, 99 )
1547
1548!
1549!-- Write buffer contents to disc immediately
1550    CALL local_flush( io )
1551
1552!
1553!-- Here the FORMATs start
1554
1555 99 FORMAT (1X,78('-'))
1556100 FORMAT (/1X,'***************************',9X,42('-')/        &
1557            1X,'* ',A,' *',9X,A/                               &
1558            1X,'***************************',9X,42('-'))
1559101 FORMAT (37X,'coupled run using MPI-',I1,': ',A/ &
1560            37X,42('-'))
1561102 FORMAT (/' Date:              ',A8,9X,'Run:       ',A20/      &
1562            ' Time:              ',A8,9X,'Run-No.:   ',I2.2/     &
1563            ' Run on host:     ',A10)
1564#if defined( __parallel )
1565103 FORMAT (' Number of PEs:',8X,I5,9X,'Processor grid (x,y): (',I3,',',I3, &
1566              ')',1X,A)
1567104 FORMAT (' Number of PEs:',8X,I5,9X,'Tasks:',I4,'   threads per task:',I4/ &
1568              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
1569105 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1570106 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
1571            37X,'because the job is running on an SMP-cluster')
1572107 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
1573#endif
1574110 FORMAT (/' Numerical Schemes:'/ &
1575             ' -----------------'/)
1576111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1577112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1578            '     Iterations (initial/other): ',I3,'/',I3,'  omega = ',F5.3)
1579113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1580                  ' or Upstream')
1581114 FORMAT (' --> Momentum advection via Upstream-Spline-Scheme')
1582115 FORMAT ('     Tendencies are smoothed via Long-Filter with factor ',F5.3) 
1583116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1584                  ' or Upstream')
1585117 FORMAT (' --> Scalar advection via Upstream-Spline-Scheme')
1586118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1587119 FORMAT (' --> Galilei-Transform applied to horizontal advection', &
1588            '     Translation velocity = ',A/ &
1589            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1590120 FORMAT (' --> Time differencing scheme: leapfrog only (no euler in case', &
1591                  ' of timestep changes)')
1592121 FORMAT (' --> Time differencing scheme: leapfrog + euler in case of', &
1593                  ' timestep changes')
1594122 FORMAT (' --> Time differencing scheme: ',A)
1595123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
1596            '     maximum damping coefficient: ',F5.3, ' 1/s')
1597124 FORMAT ('     Spline-overshoots are being suppressed')
1598125 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1599                  ' of'/                                                       &
1600            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1601126 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1602                  ' of'/                                                       &
1603            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1604127 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1605            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1606128 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1607            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1608129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1609130 FORMAT (' --> Additional prognostic equation for the total water content')
1610131 FORMAT (' --> Parameterization of condensation processes via (0%-or100%)')
1611132 FORMAT (' --> Parameterization of long-wave radiation processes via'/ &
1612            '     effective emissivity scheme')
1613133 FORMAT (' --> Precipitation parameterization via Kessler-Scheme')
1614134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1615135 FORMAT (' --> Solve perturbation pressure via multigrid method (', &
1616                  A,'-cycle)'/ &
1617            '     number of grid levels:                   ',I2/ &
1618            '     Gauss-Seidel red/black iterations:       ',I2)
1619136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1620                  I3,')')
1621137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1622            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1623                  I3,')'/ &
1624            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1625                  I3,')')
1626138 FORMAT ('     Using hybrid version for 1d-domain-decomposition')
1627139 FORMAT (' --> Loop optimization method: ',A)
1628140 FORMAT ('     maximum residual allowed:                ',E10.3)
1629141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1630142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1631                  'step')
1632143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1633                  'kinetic energy')
1634150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
1635                  'conserved'/ &
1636            '     using the ',A,' mode')
1637151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
1638152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
1639           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
1640           /'     starting from dp_level_b =', F8.3, 'm', A /)
1641153 FORMAT (' --> Large-scale vertical motion is used in the ', &
1642                  'prognostic equation for')
1643154 FORMAT ('     the potential temperature')
1644200 FORMAT (//' Run time and time step information:'/ &
1645             ' ----------------------------------'/)
1646201 FORMAT ( ' Timestep:          variable     maximum value: ',F6.3,' s', &
1647             '    CFL-factor: ',F4.2)
1648202 FORMAT ( ' Timestep:       dt = ',F6.3,' s'/)
1649203 FORMAT ( ' Start time:       ',F9.3,' s'/ &
1650             ' End time:         ',F9.3,' s')
1651204 FORMAT ( A,F9.3,' s')
1652205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
1653206 FORMAT (/' Time reached:     ',F9.3,' s'/ &
1654             ' CPU-time used:    ',F9.3,' s     per timestep:               ', &
1655               '  ',F9.3,' s'/                                                 &
1656             '                                   per second of simulated tim', &
1657               'e: ',F9.3,' s')
1658207 FORMAT ( A/' Coupling start time:',F9.3,' s')
1659250 FORMAT (//' Computational grid and domain size:'/ &
1660              ' ----------------------------------'// &
1661              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1662              ' m    dz =    ',F7.3,' m'/ &
1663              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1664              ' m  z(u) = ',F10.3,' m'/)
1665252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
1666              ' factor: ',F5.3/ &
1667            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1668254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1669            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1670255 FORMAT (' Subdomains have equal size')
1671256 FORMAT (' Subdomains at the upper edges of the virtual processor grid ', &
1672              'have smaller sizes'/                                          &
1673            ' Size of smallest subdomain:    (  ',I4,',   ',I4,',   ',I4,')')
1674260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1675             ' degrees')
1676270 FORMAT (//' Topography informations:'/ &
1677              ' -----------------------'// &
1678              1X,'Topography: ',A)
1679271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1680              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1681                ' / ',I4)
1682272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
1683              ' direction' / &
1684              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
1685              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
1686278 FORMAT (' Topography grid definition convention:'/ &
1687            ' cell edge (staggered grid points'/  &
1688            ' (u in x-direction, v in y-direction))' /)
1689279 FORMAT (' Topography grid definition convention:'/ &
1690            ' cell center (scalar grid points)' /)
1691280 FORMAT (//' Vegetation canopy (drag) model:'/ &
1692              ' ------------------------------'// &
1693              ' Canopy mode: ', A / &
1694              ' Canopy top: ',I4 / &
1695              ' Leaf drag coefficient: ',F6.2 /)
1696281 FORMAT (/ ' Scalar_exchange_coefficient: ',F6.2 / &
1697              ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /)
1698282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s')
1699283 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
1700              ' Height:              ',A,'  m'/ &
1701              ' Leaf area density:   ',A,'  m**2/m**3'/ &
1702              ' Gradient:            ',A,'  m**2/m**4'/ &
1703              ' Gridpoint:           ',A)
1704               
1705300 FORMAT (//' Boundary conditions:'/ &
1706             ' -------------------'// &
1707             '                     p                    uv             ', &
1708             '                   pt'// &
1709             ' B. bound.: ',A/ &
1710             ' T. bound.: ',A)
1711301 FORMAT (/'                     ',A// &
1712             ' B. bound.: ',A/ &
1713             ' T. bound.: ',A)
1714303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1715304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1716305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
1717               'computational u,v-level:'// &
1718             '       zp = ',F6.2,' m   z0 = ',F6.4,' m   kappa = ',F4.2/ &
1719             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1720306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
1721307 FORMAT ('       Heatflux has a random normal distribution')
1722308 FORMAT ('       Predefined surface temperature')
1723309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
1724310 FORMAT (//'    1D-Model:'// &
1725             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1726311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1727312 FORMAT ('       Predefined surface humidity')
1728313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1729314 FORMAT ('       Predefined scalar value at the surface')
1730315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
1731316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1732                    'atmosphere model')
1733317 FORMAT (//' Lateral boundaries:'/ &
1734            '       left/right:  ',A/    &
1735            '       north/south: ',A)
1736318 FORMAT (/'       outflow damping layer width: ',I3,' gridpoints with km_', &
1737                    'max =',F5.1,' m**2/s')
1738319 FORMAT ('       turbulence recycling at inflow switched on'/ &
1739            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1740            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
1741320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
1742            '                                          v: ',F9.6,' m**2/s**2')
1743325 FORMAT (//' List output:'/ &
1744             ' -----------'//  &
1745            '    1D-Profiles:'/    &
1746            '       Output every             ',F8.2,' s')
1747326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
1748            '       Averaging input every    ',F8.2,' s')
1749330 FORMAT (//' Data output:'/ &
1750             ' -----------'/)
1751331 FORMAT (/'    1D-Profiles:')
1752332 FORMAT (/'       ',A)
1753333 FORMAT ('       Output every             ',F8.2,' s',/ &
1754            '       Time averaged over       ',F8.2,' s'/ &
1755            '       Averaging input every    ',F8.2,' s')
1756334 FORMAT (/'    2D-Arrays',A,':')
1757335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1758            '       Output every             ',F8.2,' s  ',A/ &
1759            '       Cross sections at ',A1,' = ',A/ &
1760            '       scalar-coordinates:   ',A,' m'/)
1761336 FORMAT (/'    3D-Arrays',A,':')
1762337 FORMAT (/'       Arrays: ',A/ &
1763            '       Output every             ',F8.2,' s  ',A/ &
1764            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1765338 FORMAT ('       Compressed data output'/ &
1766            '       Decimal precision: ',A/)
1767339 FORMAT ('       No output during initial ',F8.2,' s')
1768340 FORMAT (/'    Time series:')
1769341 FORMAT ('       Output every             ',F8.2,' s'/)
1770342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1771            '       Output every             ',F8.2,' s  ',A/ &
1772            '       Time averaged over       ',F8.2,' s'/ &
1773            '       Averaging input every    ',F8.2,' s'/ &
1774            '       Cross sections at ',A1,' = ',A/ &
1775            '       scalar-coordinates:   ',A,' m'/)
1776343 FORMAT (/'       Arrays: ',A/ &
1777            '       Output every             ',F8.2,' s  ',A/ &
1778            '       Time averaged over       ',F8.2,' s'/ &
1779            '       Averaging input every    ',F8.2,' s'/ &
1780            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1781344 FORMAT ('       Output format: ',A/)
1782345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
1783            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
1784            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
1785            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
1786346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
1787347 FORMAT ('       Variables: ',A/ &
1788            '       Output every             ',F8.2,' s')
1789348 FORMAT ('       Variables: ',A/ &
1790            '       Output every             ',F8.2,' s'/ &
1791            '       Time averaged over       ',F8.2,' s'/ &
1792            '       Averaging input every    ',F8.2,' s')
1793349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1794            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
1795            13('       ',8(F8.2,',')/) )
1796350 FORMAT (/'       Output locations in ',A,'-direction: ', &
1797            'all gridpoints along ',A,'-direction (default).' )
1798351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1799            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
1800            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
1801#if defined( __dvrp_graphics )
1802360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
1803            '       Output every      ',F7.1,' s'/ &
1804            '       Output mode:      ',A/ &
1805            '       Host / User:      ',A,' / ',A/ &
1806            '       Directory:        ',A// &
1807            '       The sequence contains:')
1808361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
1809            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
1810362 FORMAT (/'       Slicer plane ',A/ &
1811            '       Slicer limits: [',F6.2,',',F6.2,']')
1812363 FORMAT (/'       Particles'/ &
1813            '          particle size:  ',F7.2,' m')
1814364 FORMAT ('          particle ',A,' controlled by "',A,'" with interval [', &
1815                       F6.2,',',F6.2,']')
1816365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
1817            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
1818                     ')'/ &
1819            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
1820            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
1821366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
1822367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
1823#endif
1824#if defined( __spectra )
1825370 FORMAT ('    Spectra:')
1826371 FORMAT ('       Output every ',F7.1,' s'/)
1827372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
1828            '       Directions: ', 10(A5,',')/                         &
1829            '       height levels  k = ', 20(I3,',')/                  &
1830            '                          ', 20(I3,',')/                  &
1831            '                          ', 20(I3,',')/                  &
1832            '                          ', 20(I3,',')/                  &
1833            '                          ', 19(I3,','),I3,'.'/           &
1834            '       height levels selected for standard plot:'/        &
1835            '                      k = ', 20(I3,',')/                  &
1836            '                          ', 20(I3,',')/                  &
1837            '                          ', 20(I3,',')/                  &
1838            '                          ', 20(I3,',')/                  &
1839            '                          ', 19(I3,','),I3,'.'/           &
1840            '       Time averaged over ', F7.1, ' s,' /                &
1841            '       Profiles for the time averaging are taken every ', &
1842                    F6.1,' s')
1843#endif
1844400 FORMAT (//' Physical quantities:'/ &
1845              ' -------------------'/)
1846410 FORMAT ('    Angular velocity    :   omega = ',E9.3,' rad/s'/  &
1847            '    Geograph. latitude  :   phi   = ',F4.1,' degr'/   &
1848            '    Coriolis parameter  :   f     = ',F9.6,' 1/s'/    &
1849            '                            f*    = ',F9.6,' 1/s')
1850411 FORMAT (/'    Gravity             :   g     = ',F4.1,' m/s**2')
1851412 FORMAT (/'    Reference density in buoyancy terms: ',F8.3,' kg/m**3')
1852413 FORMAT (/'    Reference temperature in buoyancy terms: ',F8.4,' K')
1853415 FORMAT (/'    Cloud physics parameters:'/ &
1854             '    ------------------------'/)
1855416 FORMAT ('        Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
1856            '        Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
1857            '        Density of air     :   rho_0 = ',F5.3,' kg/m**3'/  &
1858            '        Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
1859            '        Vapourization heat :   L_v   = ',E8.2,' J/kg')
1860420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1861            '       Height:        ',A,'  m'/ &
1862            '       Temperature:   ',A,'  K'/ &
1863            '       Gradient:      ',A,'  K/100m'/ &
1864            '       Gridpoint:     ',A)
1865421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1866            '       Height:      ',A,'  m'/ &
1867            '       Humidity:    ',A,'  kg/kg'/ &
1868            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1869            '       Gridpoint:   ',A)
1870422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1871            '       Height:                  ',A,'  m'/ &
1872            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1873            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1874            '       Gridpoint:               ',A)
1875423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1876            '       Height:      ',A,'  m'/ &
1877            '       ug:          ',A,'  m/s'/ &
1878            '       Gradient:    ',A,'  1/100s'/ &
1879            '       Gridpoint:   ',A)
1880424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1881            '       Height:      ',A,'  m'/ &
1882            '       vg:          ',A,'  m/s'/ &
1883            '       Gradient:    ',A,'  1/100s'/ &
1884            '       Gridpoint:   ',A)
1885425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
1886            '       Height:     ',A,'  m'/ &
1887            '       Salinity:   ',A,'  psu'/ &
1888            '       Gradient:   ',A,'  psu/100m'/ &
1889            '       Gridpoint:  ',A)
1890426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
1891            '       Height:      ',A,'  m'/ &
1892            '       w_subs:      ',A,'  m/s'/ &
1893            '       Gradient:    ',A,'  (m/s)/100m'/ &
1894            '       Gridpoint:   ',A)
1895450 FORMAT (//' LES / Turbulence quantities:'/ &
1896              ' ---------------------------'/)
1897451 FORMAT ('   Diffusion coefficients are constant:'/ &
1898            '   Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1899452 FORMAT ('   Mixing length is limited to the Prandtl mixing lenth.')
1900453 FORMAT ('   Mixing length is limited to ',F4.2,' * z')
1901454 FORMAT ('   TKE is not allowed to fall below ',E9.2,' (m/s)**2')
1902455 FORMAT ('   initial TKE is prescribed as ',E9.2,' (m/s)**2')
1903470 FORMAT (//' Actions during the simulation:'/ &
1904              ' -----------------------------'/)
1905471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1906            '    Disturbance amplitude           :     ',F4.2, ' m/s'/       &
1907            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1908            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
1909472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1910                 ' to i/j =',I4)
1911473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1912                 1X,F5.3, ' m**2/s**2')
1913474 FORMAT ('    Random number generator used    : ',A/)
1914475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1915                 'respectively, if'/ &
1916            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1917                 ' 3D-simulation'/)
1918476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1919                 'respectively, if the'/ &
1920            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1921                 ' the 3D-simulation'/)
1922477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1923                 'respectively, if the'/ &
1924            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1925                 ' the 3D-simulation'/)
1926480 FORMAT ('    Particles:'/ &
1927            '    ---------'// &
1928            '       Particle advection is active (switched on at t = ', F7.1, &
1929                    ' s)'/ &
1930            '       Start of new particle generations every  ',F6.1,' s'/ &
1931            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
1932            '                            bottom:     ', A, ' top:         ', A/&
1933            '       Maximum particle age:                 ',F9.1,' s'/ &
1934            '       Advection stopped at t = ',F9.1,' s'/ &
1935            '       Particles are sorted every ',F9.1,' s'/)
1936481 FORMAT ('       Particles have random start positions'/)
1937482 FORMAT ('          Particles are advected only horizontally'/)
1938483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
1939484 FORMAT ('            Number of tails of the total domain: ',I10/ &
1940            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
1941            '            Maximum age of the end of the tail:  ',F8.2,' s')
1942485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
1943486 FORMAT ('       Particle statistics are written on file'/)
1944487 FORMAT ('       Number of particle groups: ',I2/)
1945488 FORMAT ('       SGS velocity components are used for particle advection'/ &
1946            '          minimum timestep for advection: ', F7.5/)
1947489 FORMAT ('       Number of particles simultaneously released at each ', &
1948                    'point: ', I5/)
1949490 FORMAT ('       Particle group ',I2,':'/ &
1950            '          Particle radius: ',E10.3, 'm')
1951491 FORMAT ('          Particle inertia is activated'/ &
1952            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
1953492 FORMAT ('          Particles are advected only passively (no inertia)'/)
1954493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
1955            '                                         y:',F8.1,' - ',F8.1,' m'/&
1956            '                                         z:',F8.1,' - ',F8.1,' m'/&
1957            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
1958                       ' m  dz = ',F8.1,' m'/)
1959494 FORMAT ('       Output of particle time series in NetCDF format every ', &
1960                    F8.2,' s'/)
1961495 FORMAT ('       Number of particles in total domain: ',I10/)
1962500 FORMAT (//' 1D-Model parameters:'/                           &
1963              ' -------------------'//                           &
1964            '    Simulation time:                   ',F8.1,' s'/ &
1965            '    Run-controll output every:         ',F8.1,' s'/ &
1966            '    Vertical profile output every:     ',F8.1,' s'/ &
1967            '    Mixing length calculation:         ',A/         &
1968            '    Dissipation calculation:           ',A/)
1969502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
1970
1971
1972 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.