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

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

last commit documented

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