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

Last change on this file since 481 was 449, checked in by raasch, 15 years ago

branch revision comments from Marcus (rev 410) replaced by normal revision comments

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