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

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

last commit documented

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