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

Last change on this file since 668 was 668, checked in by suehring, 14 years ago

last commit documented

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