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

Last change on this file since 759 was 759, checked in by raasch, 13 years ago

New:
---

The number of parallel I/O operations can be limited with new mrun-option -w.
(advec_particles, data_output_2d, data_output_3d, header, init_grid, init_pegrid, init_3d_model, modules, palm, parin, write_3d_binary)

Changed:


mrun option -T is obligatory

Errors:


Bugfix: No zero assignments to volume_flow_initial and volume_flow_area in
case of normal restart runs. (init_3d_model)

initialization of u_0, v_0. This is just to avoid access of uninitialized
memory in exchange_horiz_2d, which causes respective error messages
when the Intel thread checker (inspector) is used. (production_e)

Bugfix for ts limitation (prandtl_fluxes)

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