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

Last change on this file since 580 was 580, checked in by heinze, 14 years ago

Renaming of ws_vertical_gradient to subs_vertical_gradient, ws_vertical_gradient_level to subs_vertical_gradient_level and ws_vertical_gradient_level_ind to subs_vertical_gradient_level_i

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