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

Last change on this file since 243 was 241, checked in by letzel, 15 years ago
  • Option to predefine a target bulk velocity for conserve_volume_flow
  • Property svn:keywords set to Id
File size: 65.5 KB
RevLine 
[1]1 SUBROUTINE header
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
[237]6! Output of cluster_size
[240]7! +canyon_height, canyon_width_x, canyon_width_y, canyon_wall_left,
[241]8! canyon_wall_south, conserve_volume_flow_mode, dp_external, dp_level_b,
9! dp_smooth, dpdxy, u_bulk, v_bulk
[1]10!
11! Former revisions:
12! -----------------
[3]13! $Id: header.f90 241 2009-02-19 16:08:56Z letzel $
[39]14!
[226]15! 206 2008-10-13 14:59:11Z raasch
16! Bugfix: error in zu index in case of section_xy = -1
17!
[200]18! 198 2008-09-17 08:55:28Z raasch
19! Format adjustments allowing output of larger revision numbers
20!
[198]21! 197 2008-09-16 15:29:03Z raasch
22! allow 100 spectra levels instead of 10 for consistency with
23! define_netcdf_header,
24! bugfix in the output of the characteristic levels of potential temperature,
25! geostrophic wind, scalar concentration, humidity and leaf area density,
26! output of turbulence recycling informations
27!
[139]28! 138 2007-11-28 10:03:58Z letzel
29! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
30! Allow two instead of one digit to specify isosurface and slicer variables.
31! Output of sorting frequency of particles
32!
[110]33! 108 2007-08-24 15:10:38Z letzel
34! Output of informations for coupled model runs (boundary conditions etc.)
35! + output of momentumfluxes at the top boundary
36! Rayleigh damping for ocean, e_init
37!
[98]38! 97 2007-06-21 08:23:15Z raasch
39! Adjustments for the ocean version.
40! use_pt_reference renamed use_reference
41!
[90]42! 87 2007-05-22 15:46:47Z raasch
43! Bugfix: output of use_upstream_for_tke
44!
[83]45! 82 2007-04-16 15:40:52Z raasch
46! Preprocessor strings for different linux clusters changed to "lc",
47! routine local_flush is used for buffer flushing
48!
[77]49! 76 2007-03-29 00:58:32Z raasch
50! Output of netcdf_64bit_3d, particles-package is now part of the default code,
51! output of the loop optimization method, moisture renamed humidity,
52! output of subversion revision number
53!
[39]54! 19 2007-02-23 04:53:48Z raasch
55! Output of scalar flux applied at top boundary
56!
[3]57! RCS Log replace by Id keyword, revision history cleaned up
58!
[1]59! Revision 1.63  2006/08/22 13:53:13  raasch
60! Output of dz_max
61!
62! Revision 1.1  1997/08/11 06:17:20  raasch
63! Initial revision
64!
65!
66! Description:
67! ------------
68! Writing a header with all important informations about the actual run.
69! This subroutine is called three times, two times at the beginning
70! (writing information on files RUN_CONTROL and HEADER) and one time at the
71! end of the run, then writing additional information about CPU-usage on file
72! header.
73!------------------------------------------------------------------------------!
74
75    USE arrays_3d
76    USE control_parameters
77    USE cloud_parameters
78    USE cpulog
79    USE dvrp_variables
80    USE grid_variables
81    USE indices
82    USE model_1d
83    USE particle_attributes
84    USE pegrid
85    USE spectrum
86
87    IMPLICIT NONE
88
89    CHARACTER (LEN=1)  ::  prec
90    CHARACTER (LEN=2)  ::  do2d_mode
91    CHARACTER (LEN=5)  ::  section_chr
92    CHARACTER (LEN=9)  ::  time_to_string
93    CHARACTER (LEN=10) ::  coor_chr, host_chr
94    CHARACTER (LEN=16) ::  begin_chr
[200]95    CHARACTER (LEN=23) ::  ver_rev
[1]96    CHARACTER (LEN=40) ::  output_format
[167]97    CHARACTER (LEN=70) ::  char1, char2, dopr_chr, &
[1]98                           do2d_xy, do2d_xz, do2d_yz, do3d_chr, &
[167]99                           run_classification
100    CHARACTER (LEN=86) ::  coordinates, gradients, learde, slices,  &
101                           temperatures, ugcomponent, vgcomponent
[1]102    CHARACTER (LEN=85) ::  roben, runten
103
[240]104    INTEGER ::  av, bh, blx, bly, bxl, bxr, byn, bys, ch, cwx, cwy, cxl, cxr, &
105         cyn, cys, i, ihost, io, j, l, ll
[1]106    REAL    ::  cpuseconds_per_simulated_second
107
108!
109!-- Open the output file. At the end of the simulation, output is directed
110!-- to unit 19.
111    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
112         .NOT. simulated_time_at_begin /= simulated_time )  THEN
113       io = 15   !  header output on file RUN_CONTROL
114    ELSE
115       io = 19   !  header output on file HEADER
116    ENDIF
117    CALL check_open( io )
118
119!
120!-- At the end of the run, output file (HEADER) will be rewritten with
121!-- new informations
122    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
123
124!
125!-- Determine kind of model run
126    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
127       run_classification = '3D - restart run'
[147]128    ELSEIF ( TRIM( initializing_actions ) == 'read_data_for_recycling' )  THEN
129       run_classification = '3D - run using 3D - prerun data'
130    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
131       run_classification = '3D - run without 1D - prerun'
[197]132    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
[147]133       run_classification = '3D - run with 1D - prerun'
[197]134    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
135       run_classification = '3D - run initialized by user'
[1]136    ELSE
[147]137       PRINT*,'+++ header:  unknown action(s): ',initializing_actions
[1]138    ENDIF
[97]139    IF ( ocean )  THEN
140       run_classification = 'ocean - ' // run_classification
141    ELSE
142       run_classification = 'atmosphere - ' // run_classification
143    ENDIF
[1]144
145!
146!-- Run-identification, date, time, host
147    host_chr = host(1:10)
[75]148    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
[102]149    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
150    IF ( coupling_mode /= 'uncoupled' )  WRITE ( io, 101 )  coupling_mode
151    WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr, &
152                       ADJUSTR( host_chr )
[1]153#if defined( __parallel )
154    IF ( npex == -1  .AND.  pdims(2) /= 1 )  THEN
155       char1 = 'calculated'
156    ELSEIF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.  &
157               host(1:2) == 'lc' )  .AND.                          &
158             npex == -1  .AND.  pdims(2) == 1 )  THEN
159       char1 = 'forced'
160    ELSE
161       char1 = 'predefined'
162    ENDIF
163    IF ( threads_per_task == 1 )  THEN
[102]164       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
[1]165    ELSE
[102]166       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
[1]167                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
168    ENDIF
169    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
170           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
171         npex == -1  .AND.  pdims(2) == 1 )                      &
172    THEN
[102]173       WRITE ( io, 106 )
[1]174    ELSEIF ( pdims(2) == 1 )  THEN
[102]175       WRITE ( io, 107 )  'x'
[1]176    ELSEIF ( pdims(1) == 1 )  THEN
[102]177       WRITE ( io, 107 )  'y'
[1]178    ENDIF
[102]179    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
[1]180#endif
181    WRITE ( io, 99 )
182
183!
184!-- Numerical schemes
185    WRITE ( io, 110 )
186    IF ( psolver(1:7) == 'poisfft' )  THEN
187       WRITE ( io, 111 )  TRIM( fft_method )
188       IF ( psolver == 'poisfft_hybrid' )  WRITE ( io, 138 )
189    ELSEIF ( psolver == 'sor' )  THEN
190       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
191    ELSEIF ( psolver == 'multigrid' )  THEN
192       WRITE ( io, 135 )  cycle_mg, maximum_grid_level, ngsrb
193       IF ( mg_cycles == -1 )  THEN
194          WRITE ( io, 140 )  residual_limit
195       ELSE
196          WRITE ( io, 141 )  mg_cycles
197       ENDIF
198       IF ( mg_switch_to_pe0_level == 0 )  THEN
199          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
200                             nzt_mg(1)
[197]201       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
[1]202          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
203                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
204                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
205                             nzt_mg(mg_switch_to_pe0_level),    &
206                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
207                             nzt_mg(1)
208       ENDIF
209    ENDIF
210    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
211    THEN
212       WRITE ( io, 142 )
213    ENDIF
214
215    IF ( momentum_advec == 'pw-scheme' )  THEN
216       WRITE ( io, 113 )
217    ELSE
218       WRITE ( io, 114 )
219       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
220       IF ( overshoot_limit_u /= 0.0  .OR.  overshoot_limit_v /= 0.0  .OR. &
221            overshoot_limit_w /= 0.0 )  THEN
222          WRITE ( io, 127 )  overshoot_limit_u, overshoot_limit_v, &
223                             overshoot_limit_w
224       ENDIF
225       IF ( ups_limit_u /= 0.0  .OR.  ups_limit_v /= 0.0  .OR. &
226            ups_limit_w /= 0.0 )                               &
227       THEN
228          WRITE ( io, 125 )  ups_limit_u, ups_limit_v, ups_limit_w
229       ENDIF
230       IF ( long_filter_factor /= 0.0 )  WRITE ( io, 115 )  long_filter_factor
231    ENDIF
232    IF ( scalar_advec == 'pw-scheme' )  THEN
233       WRITE ( io, 116 )
234    ELSEIF ( scalar_advec == 'ups-scheme' )  THEN
235       WRITE ( io, 117 )
236       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
237       IF ( overshoot_limit_e /= 0.0  .OR.  overshoot_limit_pt /= 0.0 )  THEN
238          WRITE ( io, 128 )  overshoot_limit_e, overshoot_limit_pt
239       ENDIF
240       IF ( ups_limit_e /= 0.0  .OR.  ups_limit_pt /= 0.0 )  THEN
241          WRITE ( io, 126 )  ups_limit_e, ups_limit_pt
242       ENDIF
243    ELSE
244       WRITE ( io, 118 )
245    ENDIF
[63]246
247    WRITE ( io, 139 )  TRIM( loop_optimization )
248
[1]249    IF ( galilei_transformation )  THEN
250       IF ( use_ug_for_galilei_tr )  THEN
251          char1 = 'geostrophic wind'
252       ELSE
253          char1 = 'mean wind in model domain'
254       ENDIF
255       IF ( simulated_time_at_begin == simulated_time )  THEN
256          char2 = 'at the start of the run'
257       ELSE
258          char2 = 'at the end of the run'
259       ENDIF
260       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ), &
261                          advected_distance_x/1000.0, advected_distance_y/1000.0
262    ENDIF
263    IF ( timestep_scheme == 'leapfrog' )  THEN
264       WRITE ( io, 120 )
265    ELSEIF ( timestep_scheme == 'leapfrog+euler' )  THEN
266       WRITE ( io, 121 )
267    ELSE
268       WRITE ( io, 122 )  timestep_scheme
269    ENDIF
[87]270    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
[1]271    IF ( rayleigh_damping_factor /= 0.0 )  THEN
[108]272       IF ( .NOT. ocean )  THEN
273          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
274               rayleigh_damping_factor
275       ELSE
276          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
277               rayleigh_damping_factor
278       ENDIF
[1]279    ENDIF
[75]280    IF ( humidity )  THEN
[1]281       IF ( .NOT. cloud_physics )  THEN
282          WRITE ( io, 129 )
283       ELSE
284          WRITE ( io, 130 )
285          WRITE ( io, 131 )
286          IF ( radiation )      WRITE ( io, 132 )
287          IF ( precipitation )  WRITE ( io, 133 )
288       ENDIF
289    ENDIF
290    IF ( passive_scalar )  WRITE ( io, 134 )
[240]291    IF ( conserve_volume_flow )  THEN
[241]292       WRITE ( io, 150 )  conserve_volume_flow_mode
293       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
294          WRITE ( io, 151 )  u_bulk, v_bulk
295       ENDIF
[240]296    ELSEIF ( dp_external )  THEN
297       IF ( dp_smooth )  THEN
[241]298          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
[240]299       ELSE
[241]300          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
[240]301       ENDIF
302    ENDIF
[1]303    WRITE ( io, 99 )
304
305!
306!-- Runtime and timestep informations
307    WRITE ( io, 200 )
308    IF ( .NOT. dt_fixed )  THEN
309       WRITE ( io, 201 )  dt_max, cfl_factor
310    ELSE
311       WRITE ( io, 202 )  dt
312    ENDIF
313    WRITE ( io, 203 )  simulated_time_at_begin, end_time
314
315    IF ( time_restart /= 9999999.9  .AND. &
316         simulated_time_at_begin == simulated_time )  THEN
317       IF ( dt_restart == 9999999.9 )  THEN
318          WRITE ( io, 204 )  ' Restart at:       ',time_restart
319       ELSE
320          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
321       ENDIF
322    ENDIF
323
324    IF ( simulated_time_at_begin /= simulated_time )  THEN
325       i = MAX ( log_point_s(10)%counts, 1 )
326       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0 )  THEN
327          cpuseconds_per_simulated_second = 0.0
328       ELSE
329          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
330                                            ( simulated_time -    &
331                                              simulated_time_at_begin )
332       ENDIF
333       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum, &
334                          log_point_s(10)%sum / REAL( i ),     &
335                          cpuseconds_per_simulated_second
336       IF ( time_restart /= 9999999.9  .AND.  time_restart < end_time )  THEN
337          IF ( dt_restart == 9999999.9 )  THEN
338             WRITE ( io, 204 )  ' Next restart at:  ',time_restart
339          ELSE
340             WRITE ( io, 205 )  ' Next restart at:  ',time_restart, dt_restart
341          ENDIF
342       ENDIF
343    ENDIF
344
345!
346!-- Computational grid
[94]347    IF ( .NOT. ocean )  THEN
348       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
349       IF ( dz_stretch_level_index < nzt+1 )  THEN
350          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
351                             dz_stretch_factor, dz_max
352       ENDIF
353    ELSE
354       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
355       IF ( dz_stretch_level_index > 0 )  THEN
356          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
357                             dz_stretch_factor, dz_max
358       ENDIF
[1]359    ENDIF
360    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
361                       MIN( nnz+2, nzt+2 )
[76]362    IF ( numprocs > 1 )  THEN
363       IF ( nxa == nx  .AND.  nya == ny  .AND.  nza == nz )  THEN
364          WRITE ( io, 255 )
365       ELSE
366          WRITE ( io, 256 )  nnx-(nxa-nx), nny-(nya-ny), nzt+2
367       ENDIF
[1]368    ENDIF
369    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
370
371!
372!-- Topography
373    WRITE ( io, 270 )  topography
374    SELECT CASE ( TRIM( topography ) )
375
376       CASE ( 'flat' )
377          ! no actions necessary
378
379       CASE ( 'single_building' )
380          blx = INT( building_length_x / dx )
381          bly = INT( building_length_y / dy )
382          bh  = INT( building_height / dz )
383
384          IF ( building_wall_left == 9999999.9 )  THEN
385             building_wall_left = ( nx + 1 - blx ) / 2 * dx
386          ENDIF
387          bxl = INT ( building_wall_left / dx + 0.5 )
388          bxr = bxl + blx
389
390          IF ( building_wall_south == 9999999.9 )  THEN
391             building_wall_south = ( ny + 1 - bly ) / 2 * dy
392          ENDIF
393          bys = INT ( building_wall_south / dy + 0.5 )
394          byn = bys + bly
395
396          WRITE ( io, 271 )  building_length_x, building_length_y, &
397                             building_height, bxl, bxr, bys, byn
398
[240]399       CASE ( 'single_street_canyon' )
400          ch  = NINT( canyon_height / dz )
401          IF ( canyon_width_x /= 9999999.9 )  THEN
402!
403!--          Street canyon in y direction
404             cwx = NINT( canyon_width_x / dx )
405             IF ( canyon_wall_left == 9999999.9 )  THEN
406                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
407             ENDIF
408             cxl = NINT( canyon_wall_left / dx )
409             cxr = cxl + cwx
410             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
411
412          ELSEIF ( canyon_width_y /= 9999999.9 )  THEN
413!
414!--          Street canyon in x direction
415             cwy = NINT( canyon_width_y / dy )
416             IF ( canyon_wall_south == 9999999.9 )  THEN
417                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
418             ENDIF
419             cys = NINT( canyon_wall_south / dy )
420             cyn = cys + cwy
421             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
422          ENDIF
423
[1]424    END SELECT
425
[138]426    IF ( plant_canopy ) THEN
427
428       WRITE ( io, 280 ) canopy_mode, pch_index, drag_coefficient
[153]429       IF ( passive_scalar ) THEN
430          WRITE ( io, 281 ) scalar_exchange_coefficient,   &
431                            leaf_surface_concentration
432       ENDIF
[138]433
[1]434!
[153]435!--    Heat flux at the top of vegetation
436       WRITE ( io, 282 ) cthf
437
438!
[138]439!--    Leaf area density profile
440!--    Building output strings, starting with surface value
441       WRITE ( learde, '(F6.2)' )  lad_surface
442       gradients = '------'
443       slices = '     0'
444       coordinates = '   0.0'
445       i = 1
446       DO  WHILE ( lad_vertical_gradient_level_ind(i) /= -9999 )
447
448          WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
449          learde = TRIM( learde ) // ' ' // TRIM( coor_chr )
450
451          WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
452          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
453
454          WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
455          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
456
457          WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
458          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
459
460          i = i + 1
461       ENDDO
462
[153]463       WRITE ( io, 283 )  TRIM( coordinates ), TRIM( learde ), &
[138]464                          TRIM( gradients ), TRIM( slices )
465
466    ENDIF
467
468!
[1]469!-- Boundary conditions
470    IF ( ibc_p_b == 0 )  THEN
471       runten = 'p(0)     = 0      |'
472    ELSEIF ( ibc_p_b == 1 )  THEN
473       runten = 'p(0)     = p(1)   |'
474    ELSE
475       runten = 'p(0)     = p(1) +R|'
476    ENDIF
477    IF ( ibc_p_t == 0 )  THEN
478       roben  = 'p(nzt+1) = 0      |'
479    ELSE
480       roben  = 'p(nzt+1) = p(nzt) |'
481    ENDIF
482
483    IF ( ibc_uv_b == 0 )  THEN
484       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
485    ELSE
486       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
487    ENDIF
[132]488    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
489       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
490    ELSEIF ( ibc_uv_t == 0 )  THEN
[1]491       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
492    ELSE
493       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
494    ENDIF
495
496    IF ( ibc_pt_b == 0 )  THEN
497       runten = TRIM( runten ) // ' pt(0)   = pt_surface'
[102]498    ELSEIF ( ibc_pt_b == 1 )  THEN
[1]499       runten = TRIM( runten ) // ' pt(0)   = pt(1)'
[102]500    ELSEIF ( ibc_pt_b == 2 )  THEN
501       runten = TRIM( runten ) // ' pt(0) = from coupled model'
[1]502    ENDIF
503    IF ( ibc_pt_t == 0 )  THEN
[19]504       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
505    ELSEIF( ibc_pt_t == 1 )  THEN
506       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
507    ELSEIF( ibc_pt_t == 2 )  THEN
508       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
[1]509    ENDIF
510
511    WRITE ( io, 300 )  runten, roben
512
513    IF ( .NOT. constant_diffusion )  THEN
514       IF ( ibc_e_b == 1 )  THEN
515          runten = 'e(0)     = e(1)'
516       ELSE
517          runten = 'e(0)     = e(1) = (u*/0.1)**2'
518       ENDIF
519       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
520
[97]521       WRITE ( io, 301 )  'e', runten, roben       
[1]522
523    ENDIF
524
[97]525    IF ( ocean )  THEN
526       runten = 'sa(0)    = sa(1)'
527       IF ( ibc_sa_t == 0 )  THEN
528          roben =  'sa(nzt+1) = sa_surface'
[1]529       ELSE
[97]530          roben =  'sa(nzt+1) = sa(nzt)'
[1]531       ENDIF
[97]532       WRITE ( io, 301 ) 'sa', runten, roben
533    ENDIF
[1]534
[97]535    IF ( humidity )  THEN
536       IF ( ibc_q_b == 0 )  THEN
537          runten = 'q(0)     = q_surface'
538       ELSE
539          runten = 'q(0)     = q(1)'
540       ENDIF
541       IF ( ibc_q_t == 0 )  THEN
542          roben =  'q(nzt)   = q_top'
543       ELSE
544          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
545       ENDIF
546       WRITE ( io, 301 ) 'q', runten, roben
547    ENDIF
[1]548
[97]549    IF ( passive_scalar )  THEN
550       IF ( ibc_q_b == 0 )  THEN
551          runten = 's(0)     = s_surface'
552       ELSE
553          runten = 's(0)     = s(1)'
554       ENDIF
555       IF ( ibc_q_t == 0 )  THEN
556          roben =  's(nzt)   = s_top'
557       ELSE
558          roben =  's(nzt)   = s(nzt-1) + ds/dz'
559       ENDIF
560       WRITE ( io, 301 ) 's', runten, roben
[1]561    ENDIF
562
563    IF ( use_surface_fluxes )  THEN
564       WRITE ( io, 303 )
565       IF ( constant_heatflux )  THEN
566          WRITE ( io, 306 )  surface_heatflux
567          IF ( random_heatflux )  WRITE ( io, 307 )
568       ENDIF
[75]569       IF ( humidity  .AND.  constant_waterflux )  THEN
[1]570          WRITE ( io, 311 ) surface_waterflux
571       ENDIF
572       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
573          WRITE ( io, 313 ) surface_waterflux
574       ENDIF
575    ENDIF
576
[19]577    IF ( use_top_fluxes )  THEN
578       WRITE ( io, 304 )
[102]579       IF ( coupling_mode == 'uncoupled' )  THEN
[151]580          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
[102]581          IF ( constant_top_heatflux )  THEN
582             WRITE ( io, 306 )  top_heatflux
583          ENDIF
584       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
585          WRITE ( io, 316 )
[19]586       ENDIF
[97]587       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
588          WRITE ( io, 309 )  top_salinityflux
589       ENDIF
[75]590       IF ( humidity  .OR.  passive_scalar )  THEN
[19]591          WRITE ( io, 315 )
592       ENDIF
593    ENDIF
594
[1]595    IF ( prandtl_layer )  THEN
[94]596       WRITE ( io, 305 )  0.5 * (zu(1)-zu(0)), roughness_length, kappa, &
597                          rif_min, rif_max
[1]598       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
[75]599       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
[1]600          WRITE ( io, 312 )
601       ENDIF
602       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
603          WRITE ( io, 314 )
604       ENDIF
605    ELSE
606       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
607          WRITE ( io, 310 )  rif_min, rif_max
608       ENDIF
609    ENDIF
610
611    WRITE ( io, 317 )  bc_lr, bc_ns
612    IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
613       WRITE ( io, 318 )  outflow_damping_width, km_damp_max
[151]614       IF ( turbulent_inflow )  THEN
615          WRITE ( io, 319 )  recycling_width, recycling_plane, &
616                             inflow_damping_height, inflow_damping_width
617       ENDIF
[1]618    ENDIF
619
620!
621!-- Listing of 1D-profiles
[151]622    WRITE ( io, 325 )  dt_dopr_listing
[1]623    IF ( averaging_interval_pr /= 0.0 )  THEN
[151]624       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]625    ENDIF
626
627!
628!-- DATA output
629    WRITE ( io, 330 )
630    IF ( averaging_interval_pr /= 0.0 )  THEN
[151]631       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]632    ENDIF
633
634!
635!-- 1D-profiles
636    dopr_chr = 'Profile:'
637    IF ( dopr_n /= 0 )  THEN
638       WRITE ( io, 331 )
639
640       output_format = ''
641       IF ( netcdf_output )  THEN
642          IF ( netcdf_64bit )  THEN
643             output_format = 'netcdf (64 bit offset)'
644          ELSE
645             output_format = 'netcdf'
646          ENDIF
647       ENDIF
648       IF ( profil_output )  THEN
649          IF ( netcdf_output )  THEN
650             output_format = TRIM( output_format ) // ' and profil'
651          ELSE
652             output_format = 'profil'
653          ENDIF
654       ENDIF
655       WRITE ( io, 345 )  output_format
656
657       DO  i = 1, dopr_n
658          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
659          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
660             WRITE ( io, 332 )  dopr_chr
661             dopr_chr = '       :'
662          ENDIF
663       ENDDO
664
665       IF ( dopr_chr /= '' )  THEN
666          WRITE ( io, 332 )  dopr_chr
667       ENDIF
668       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
669       IF ( skip_time_dopr /= 0.0 )  WRITE ( io, 339 )  skip_time_dopr
670    ENDIF
671
672!
673!-- 2D-arrays
674    DO  av = 0, 1
675
676       i = 1
677       do2d_xy = ''
678       do2d_xz = ''
679       do2d_yz = ''
680       DO  WHILE ( do2d(av,i) /= ' ' )
681
682          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
683          do2d_mode = do2d(av,i)(l-1:l)
684
685          SELECT CASE ( do2d_mode )
686             CASE ( 'xy' )
687                ll = LEN_TRIM( do2d_xy )
688                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
689             CASE ( 'xz' )
690                ll = LEN_TRIM( do2d_xz )
691                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
692             CASE ( 'yz' )
693                ll = LEN_TRIM( do2d_yz )
694                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
695          END SELECT
696
697          i = i + 1
698
699       ENDDO
700
701       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
702              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
703              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) )  .AND. &
704            ( netcdf_output  .OR.  iso2d_output ) )  THEN
705
706          IF (  av == 0 )  THEN
707             WRITE ( io, 334 )  ''
708          ELSE
709             WRITE ( io, 334 )  '(time-averaged)'
710          ENDIF
711
712          IF ( do2d_at_begin )  THEN
713             begin_chr = 'and at the start'
714          ELSE
715             begin_chr = ''
716          ENDIF
717
718          output_format = ''
719          IF ( netcdf_output )  THEN
720             IF ( netcdf_64bit )  THEN
721                output_format = 'netcdf (64 bit offset)'
722             ELSE
723                output_format = 'netcdf'
724             ENDIF
725          ENDIF
726          IF ( iso2d_output )  THEN
727             IF ( netcdf_output )  THEN
728                output_format = TRIM( output_format ) // ' and iso2d'
729             ELSE
730                output_format = 'iso2d'
731             ENDIF
732          ENDIF
733          WRITE ( io, 345 )  output_format
734
735          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
736             i = 1
737             slices = '/'
738             coordinates = '/'
739!
740!--          Building strings with index and coordinate informations of the
741!--          slices
742             DO  WHILE ( section(i,1) /= -9999 )
743
744                WRITE (section_chr,'(I5)')  section(i,1)
745                section_chr = ADJUSTL( section_chr )
746                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
747
[206]748                IF ( section(i,1) == -1 )  THEN
749                   WRITE (coor_chr,'(F10.1)')  -1.0
750                ELSE
751                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
752                ENDIF
[1]753                coor_chr = ADJUSTL( coor_chr )
754                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
755
756                i = i + 1
757             ENDDO
758             IF ( av == 0 )  THEN
759                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
760                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
761                                   TRIM( coordinates )
762                IF ( skip_time_do2d_xy /= 0.0 )  THEN
763                   WRITE ( io, 339 )  skip_time_do2d_xy
764                ENDIF
765             ELSE
766                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
767                                   TRIM( begin_chr ), averaging_interval, &
768                                   dt_averaging_input, 'k', TRIM( slices ), &
769                                   TRIM( coordinates )
770                IF ( skip_time_data_output_av /= 0.0 )  THEN
771                   WRITE ( io, 339 )  skip_time_data_output_av
772                ENDIF
773             ENDIF
774
775          ENDIF
776
777          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
778             i = 1
779             slices = '/'
780             coordinates = '/'
781!
782!--          Building strings with index and coordinate informations of the
783!--          slices
784             DO  WHILE ( section(i,2) /= -9999 )
785
786                WRITE (section_chr,'(I5)')  section(i,2)
787                section_chr = ADJUSTL( section_chr )
788                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
789
790                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
791                coor_chr = ADJUSTL( coor_chr )
792                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
793
794                i = i + 1
795             ENDDO
796             IF ( av == 0 )  THEN
797                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
798                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
799                                   TRIM( coordinates )
800                IF ( skip_time_do2d_xz /= 0.0 )  THEN
801                   WRITE ( io, 339 )  skip_time_do2d_xz
802                ENDIF
803             ELSE
804                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
805                                   TRIM( begin_chr ), averaging_interval, &
806                                   dt_averaging_input, 'j', TRIM( slices ), &
807                                   TRIM( coordinates )
808                IF ( skip_time_data_output_av /= 0.0 )  THEN
809                   WRITE ( io, 339 )  skip_time_data_output_av
810                ENDIF
811             ENDIF
812          ENDIF
813
814          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
815             i = 1
816             slices = '/'
817             coordinates = '/'
818!
819!--          Building strings with index and coordinate informations of the
820!--          slices
821             DO  WHILE ( section(i,3) /= -9999 )
822
823                WRITE (section_chr,'(I5)')  section(i,3)
824                section_chr = ADJUSTL( section_chr )
825                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
826
827                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
828                coor_chr = ADJUSTL( coor_chr )
829                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
830
831                i = i + 1
832             ENDDO
833             IF ( av == 0 )  THEN
834                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
835                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
836                                   TRIM( coordinates )
837                IF ( skip_time_do2d_yz /= 0.0 )  THEN
838                   WRITE ( io, 339 )  skip_time_do2d_yz
839                ENDIF
840             ELSE
841                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
842                                   TRIM( begin_chr ), averaging_interval, &
843                                   dt_averaging_input, 'i', TRIM( slices ), &
844                                   TRIM( coordinates )
845                IF ( skip_time_data_output_av /= 0.0 )  THEN
846                   WRITE ( io, 339 )  skip_time_data_output_av
847                ENDIF
848             ENDIF
849          ENDIF
850
851       ENDIF
852
853    ENDDO
854
855!
856!-- 3d-arrays
857    DO  av = 0, 1
858
859       i = 1
860       do3d_chr = ''
861       DO  WHILE ( do3d(av,i) /= ' ' )
862
863          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
864          i = i + 1
865
866       ENDDO
867
868       IF ( do3d_chr /= '' )  THEN
869          IF ( av == 0 )  THEN
870             WRITE ( io, 336 )  ''
871          ELSE
872             WRITE ( io, 336 )  '(time-averaged)'
873          ENDIF
874
875          output_format = ''
876          IF ( netcdf_output )  THEN
[46]877             IF ( netcdf_64bit .AND. netcdf_64bit_3d )  THEN
[1]878                output_format = 'netcdf (64 bit offset)'
879             ELSE
880                output_format = 'netcdf'
881             ENDIF
882          ENDIF
883          IF ( avs_output )  THEN
884             IF ( netcdf_output )  THEN
885                output_format = TRIM( output_format ) // ' and avs'
886             ELSE
887                output_format = 'avs'
888             ENDIF
889          ENDIF
890          WRITE ( io, 345 )  output_format
891
892          IF ( do3d_at_begin )  THEN
893             begin_chr = 'and at the start'
894          ELSE
895             begin_chr = ''
896          ENDIF
897          IF ( av == 0 )  THEN
898             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
899                                zu(nz_do3d), nz_do3d
900          ELSE
901             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
902                                TRIM( begin_chr ), averaging_interval, &
903                                dt_averaging_input, zu(nz_do3d), nz_do3d
904          ENDIF
905
906          IF ( do3d_compress )  THEN
907             do3d_chr = ''
908             i = 1
909             DO WHILE ( do3d(av,i) /= ' ' )
910
911                SELECT CASE ( do3d(av,i) )
912                   CASE ( 'u' )
913                      j = 1
914                   CASE ( 'v' )
915                      j = 2
916                   CASE ( 'w' )
917                      j = 3
918                   CASE ( 'p' )
919                      j = 4
920                   CASE ( 'pt' )
921                      j = 5
922                END SELECT
923                WRITE ( prec, '(I1)' )  plot_3d_precision(j)%precision
924                do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // &
925                           ':' // prec // ','
926                i = i + 1
927
928             ENDDO
929             WRITE ( io, 338 )  do3d_chr
930
931          ENDIF
932
933          IF ( av == 0 )  THEN
934             IF ( skip_time_do3d /= 0.0 )  THEN
935                WRITE ( io, 339 )  skip_time_do3d
936             ENDIF
937          ELSE
938             IF ( skip_time_data_output_av /= 0.0 )  THEN
939                WRITE ( io, 339 )  skip_time_data_output_av
940             ENDIF
941          ENDIF
942
943       ENDIF
944
945    ENDDO
946
947!
948!-- Timeseries
949    IF ( dt_dots /= 9999999.9 )  THEN
950       WRITE ( io, 340 )
951
952       output_format = ''
953       IF ( netcdf_output )  THEN
954          IF ( netcdf_64bit )  THEN
955             output_format = 'netcdf (64 bit offset)'
956          ELSE
957             output_format = 'netcdf'
958          ENDIF
959       ENDIF
960       IF ( profil_output )  THEN
961          IF ( netcdf_output )  THEN
962             output_format = TRIM( output_format ) // ' and profil'
963          ELSE
964             output_format = 'profil'
965          ENDIF
966       ENDIF
967       WRITE ( io, 345 )  output_format
968       WRITE ( io, 341 )  dt_dots
969    ENDIF
970
971#if defined( __dvrp_graphics )
972!
973!-- Dvrp-output
974    IF ( dt_dvrp /= 9999999.9 )  THEN
975       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
976                          TRIM( dvrp_username ), TRIM( dvrp_directory )
977       i = 1
978       l = 0
979       DO WHILE ( mode_dvrp(i) /= ' ' )
980          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
[130]981             READ ( mode_dvrp(i), '(10X,I2)' )  j
[1]982             l = l + 1
983             IF ( do3d(0,j) /= ' ' )  THEN
984                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l)
985             ENDIF
986          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
[130]987             READ ( mode_dvrp(i), '(6X,I2)' )  j
[1]988             IF ( do2d(0,j) /= ' ' )  WRITE ( io, 362 )  TRIM( do2d(0,j) )
989          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
990             WRITE ( io, 363 )
991          ENDIF
992          i = i + 1
993       ENDDO
[237]994
995       IF ( TRIM( topography ) /= 'flat'  .AND.  cluster_size > 1 )  THEN
996          WRITE ( io, 364 )  cluster_size
997       ENDIF
998
[1]999    ENDIF
1000#endif
1001
1002#if defined( __spectra )
1003!
1004!-- Spectra output
1005    IF ( dt_dosp /= 9999999.9 ) THEN
1006       WRITE ( io, 370 )
1007
1008       output_format = ''
1009       IF ( netcdf_output )  THEN
1010          IF ( netcdf_64bit )  THEN
1011             output_format = 'netcdf (64 bit offset)'
1012          ELSE
1013             output_format = 'netcdf'
1014          ENDIF
1015       ENDIF
1016       IF ( profil_output )  THEN
1017          IF ( netcdf_output )  THEN
1018             output_format = TRIM( output_format ) // ' and profil'
1019          ELSE
1020             output_format = 'profil'
1021          ENDIF
1022       ENDIF
1023       WRITE ( io, 345 )  output_format
1024       WRITE ( io, 371 )  dt_dosp
1025       IF ( skip_time_dosp /= 0.0 )  WRITE ( io, 339 )  skip_time_dosp
1026       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
1027                          ( spectra_direction(i), i = 1,10 ),  &
[189]1028                          ( comp_spectra_level(i), i = 1,100 ), &
1029                          ( plot_spectra_level(i), i = 1,100 ), &
[1]1030                          averaging_interval_sp, dt_averaging_input_pr
1031    ENDIF
1032#endif
1033
1034    WRITE ( io, 99 )
1035
1036!
1037!-- Physical quantities
1038    WRITE ( io, 400 )
1039
1040!
1041!-- Geostrophic parameters
1042    WRITE ( io, 410 )  omega, phi, f, fs
1043
1044!
1045!-- Other quantities
1046    WRITE ( io, 411 )  g
[97]1047    IF ( use_reference )  THEN
1048       IF ( ocean )  THEN
1049          WRITE ( io, 412 )  prho_reference
1050       ELSE
1051          WRITE ( io, 413 )  pt_reference
1052       ENDIF
1053    ENDIF
[1]1054
1055!
1056!-- Cloud physics parameters
1057    IF ( cloud_physics ) THEN
[57]1058       WRITE ( io, 415 )
1059       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
[1]1060    ENDIF
1061
1062!-- Profile of the geostrophic wind (component ug)
1063!-- Building output strings
1064    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
1065    gradients = '------'
1066    slices = '     0'
1067    coordinates = '   0.0'
1068    i = 1
1069    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
1070     
[167]1071       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
[1]1072       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
1073
[167]1074       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
[1]1075       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1076
[167]1077       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
[1]1078       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1079
[167]1080       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
[1]1081       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1082
1083       i = i + 1
1084    ENDDO
1085
1086    WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
1087                       TRIM( gradients ), TRIM( slices )
1088
1089!-- Profile of the geostrophic wind (component vg)
1090!-- Building output strings
1091    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
1092    gradients = '------'
1093    slices = '     0'
1094    coordinates = '   0.0'
1095    i = 1
1096    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
1097
[167]1098       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
[1]1099       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
1100
[167]1101       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
[1]1102       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1103
[167]1104       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
[1]1105       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1106
[167]1107       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
[1]1108       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1109
1110       i = i + 1 
1111    ENDDO
1112
1113    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
1114                       TRIM( gradients ), TRIM( slices )
1115
1116!
1117!-- Initial temperature profile
1118!-- Building output strings, starting with surface temperature
1119    WRITE ( temperatures, '(F6.2)' )  pt_surface
1120    gradients = '------'
1121    slices = '     0'
1122    coordinates = '   0.0'
1123    i = 1
1124    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1125
[94]1126       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1127       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
[1]1128
[94]1129       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1130       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
[1]1131
[94]1132       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1133       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
[1]1134
[94]1135       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1136       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
[1]1137
1138       i = i + 1
1139    ENDDO
1140
1141    WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1142                       TRIM( gradients ), TRIM( slices )
1143
1144!
1145!-- Initial humidity profile
1146!-- Building output strings, starting with surface humidity
[75]1147    IF ( humidity  .OR.  passive_scalar )  THEN
[1]1148       WRITE ( temperatures, '(E8.1)' )  q_surface
1149       gradients = '--------'
1150       slices = '       0'
1151       coordinates = '     0.0'
1152       i = 1
1153       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1154         
1155          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1156          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1157
1158          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1159          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1160         
1161          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1162          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1163         
1164          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1165          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1166
1167          i = i + 1
1168       ENDDO
1169
[75]1170       IF ( humidity )  THEN
[1]1171          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1172                             TRIM( gradients ), TRIM( slices )
1173       ELSE
1174          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1175                             TRIM( gradients ), TRIM( slices )
1176       ENDIF
1177    ENDIF
1178
1179!
[97]1180!-- Initial salinity profile
1181!-- Building output strings, starting with surface salinity
1182    IF ( ocean )  THEN
1183       WRITE ( temperatures, '(F6.2)' )  sa_surface
1184       gradients = '------'
1185       slices = '     0'
1186       coordinates = '   0.0'
1187       i = 1
1188       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1189
1190          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1191          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1192
1193          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1194          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1195
1196          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1197          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1198
1199          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1200          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1201
1202          i = i + 1
1203       ENDDO
1204
1205       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1206                          TRIM( gradients ), TRIM( slices )
1207    ENDIF
1208
1209!
[1]1210!-- LES / turbulence parameters
1211    WRITE ( io, 450 )
1212
1213!--
1214! ... LES-constants used must still be added here
1215!--
1216    IF ( constant_diffusion )  THEN
1217       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1218                          prandtl_number
1219    ENDIF
1220    IF ( .NOT. constant_diffusion)  THEN
[108]1221       IF ( e_init > 0.0 )  WRITE ( io, 455 )  e_init
[1]1222       IF ( e_min > 0.0 )  WRITE ( io, 454 )  e_min
1223       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1224       IF ( adjust_mixing_length  .AND.  prandtl_layer )  WRITE ( io, 452 )
1225    ENDIF
1226
1227!
1228!-- Special actions during the run
1229    WRITE ( io, 470 )
1230    IF ( create_disturbances )  THEN
1231       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1232                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1233                          zu(disturbance_level_ind_t), disturbance_level_ind_t
1234       IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
1235          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1236       ELSE
1237          WRITE ( io, 473 )  disturbance_energy_limit
1238       ENDIF
1239       WRITE ( io, 474 )  TRIM( random_generator )
1240    ENDIF
1241    IF ( pt_surface_initial_change /= 0.0 )  THEN
1242       WRITE ( io, 475 )  pt_surface_initial_change
1243    ENDIF
[75]1244    IF ( humidity  .AND.  q_surface_initial_change /= 0.0 )  THEN
[1]1245       WRITE ( io, 476 )  q_surface_initial_change       
1246    ENDIF
1247    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0 )  THEN
1248       WRITE ( io, 477 )  q_surface_initial_change       
1249    ENDIF
1250
[60]1251    IF ( particle_advection )  THEN
[1]1252!
[60]1253!--    Particle attributes
1254       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1255                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
[117]1256                          end_time_prel, dt_sort_particles
[60]1257       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1258       IF ( random_start_position )  WRITE ( io, 481 )
1259       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1260       WRITE ( io, 495 )  total_number_of_particles
1261       IF ( .NOT. vertical_particle_advection )  WRITE ( io, 482 )
1262       IF ( maximum_number_of_tailpoints /= 0 )  THEN
1263          WRITE ( io, 483 )  maximum_number_of_tailpoints
1264          IF ( minimum_tailpoint_distance /= 0 )  THEN
1265             WRITE ( io, 484 )  total_number_of_tails,      &
1266                                minimum_tailpoint_distance, &
1267                                maximum_tailpoint_age
1268          ENDIF
[1]1269       ENDIF
[60]1270       IF ( dt_write_particle_data /= 9999999.9 )  THEN
1271          WRITE ( io, 485 )  dt_write_particle_data
1272          output_format = ''
1273          IF ( netcdf_output )  THEN
1274             IF ( netcdf_64bit )  THEN
1275                output_format = 'netcdf (64 bit offset) and binary'
1276             ELSE
1277                output_format = 'netcdf and binary'
1278             ENDIF
[1]1279          ELSE
[60]1280             output_format = 'binary'
[1]1281          ENDIF
[60]1282          WRITE ( io, 345 )  output_format
[1]1283       ENDIF
[60]1284       IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1285       IF ( write_particle_statistics )  WRITE ( io, 486 )
[1]1286
[60]1287       WRITE ( io, 487 )  number_of_particle_groups
[1]1288
[60]1289       DO  i = 1, number_of_particle_groups
1290          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9 )  THEN
1291             WRITE ( io, 490 )  i, 0.0
1292             WRITE ( io, 492 )
[1]1293          ELSE
[60]1294             WRITE ( io, 490 )  i, radius(i)
1295             IF ( density_ratio(i) /= 0.0 )  THEN
1296                WRITE ( io, 491 )  density_ratio(i)
1297             ELSE
1298                WRITE ( io, 492 )
1299             ENDIF
[1]1300          ENDIF
[60]1301          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1302                             pdx(i), pdy(i), pdz(i)
1303       ENDDO
[1]1304
[60]1305    ENDIF
[1]1306
[60]1307
[1]1308!
1309!-- Parameters of 1D-model
1310    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1311       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1312                          mixing_length_1d, dissipation_1d
1313       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1314          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1315       ENDIF
1316    ENDIF
1317
1318!
1319!-- User-defined informations
1320    CALL user_header( io )
1321
1322    WRITE ( io, 99 )
1323
1324!
1325!-- Write buffer contents to disc immediately
[82]1326    CALL local_flush( io )
[1]1327
1328!
1329!-- Here the FORMATs start
1330
1331 99 FORMAT (1X,78('-'))
[200]1332100 FORMAT (/1X,'***************************',9X,42('-')/        &
1333            1X,'* ',A,' *',9X,A/                               &
1334            1X,'***************************',9X,42('-'))
[102]1335101 FORMAT (37X,'coupled run: ',A/ &
1336            37X,42('-'))
[200]1337102 FORMAT (/' Date:              ',A8,9X,'Run:       ',A20/      &
1338            ' Time:              ',A8,9X,'Run-No.:   ',I2.2/     &
1339            ' Run on host:     ',A10)
[1]1340#if defined( __parallel )
[200]1341103 FORMAT (' Number of PEs:',8X,I5,9X,'Processor grid (x,y): (',I3,',',I3, &
[1]1342              ')',1X,A)
[200]1343104 FORMAT (' Number of PEs:',8X,I5,9X,'Tasks:',I4,'   threads per task:',I4/ &
[1]1344              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
[102]1345105 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1346106 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
[1]1347            37X,'because the job is running on an SMP-cluster')
[102]1348107 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
[1]1349#endif
1350110 FORMAT (/' Numerical Schemes:'/ &
1351             ' -----------------'/)
1352111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1353112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1354            '     Iterations (initial/other): ',I3,'/',I3,'  omega = ',F5.3)
1355113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1356                  ' or Upstream')
1357114 FORMAT (' --> Momentum advection via Upstream-Spline-Scheme')
1358115 FORMAT ('     Tendencies are smoothed via Long-Filter with factor ',F5.3) 
1359116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1360                  ' or Upstream')
1361117 FORMAT (' --> Scalar advection via Upstream-Spline-Scheme')
1362118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1363119 FORMAT (' --> Galilei-Transform applied to horizontal advection', &
1364            '     Translation velocity = ',A/ &
1365            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1366120 FORMAT (' --> Time differencing scheme: leapfrog only (no euler in case', &
1367                  ' of timestep changes)')
1368121 FORMAT (' --> Time differencing scheme: leapfrog + euler in case of', &
1369                  ' timestep changes')
1370122 FORMAT (' --> Time differencing scheme: ',A)
[108]1371123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
[1]1372            '     maximum damping coefficient: ',F5.3, ' 1/s')
1373124 FORMAT ('     Spline-overshoots are being suppressed')
1374125 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1375                  ' of'/                                                       &
1376            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1377126 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1378                  ' of'/                                                       &
1379            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1380127 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1381            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1382128 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1383            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1384129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1385130 FORMAT (' --> Additional prognostic equation for the total water content')
1386131 FORMAT (' --> Parameterization of condensation processes via (0%-or100%)')
1387132 FORMAT (' --> Parameterization of long-wave radiation processes via'/ &
1388            '     effective emissivity scheme')
1389133 FORMAT (' --> Precipitation parameterization via Kessler-Scheme')
1390134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1391135 FORMAT (' --> Solve perturbation pressure via multigrid method (', &
1392                  A,'-cycle)'/ &
1393            '     number of grid levels:                   ',I2/ &
1394            '     Gauss-Seidel red/black iterations:       ',I2)
1395136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1396                  I3,')')
1397137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1398            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1399                  I3,')'/ &
1400            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1401                  I3,')')
1402138 FORMAT ('     Using hybrid version for 1d-domain-decomposition')
[63]1403139 FORMAT (' --> Loop optimization method: ',A)
[1]1404140 FORMAT ('     maximum residual allowed:                ',E10.3)
1405141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1406142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1407                  'step')
[87]1408143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1409                  'kinetic energy')
[1]1410150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
[241]1411                  'conserved'/ &
1412            '     using the ',A,' mode')
1413151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
1414152 FORMAT (' --> External pressure gradient directly prescribed by the user:'/, &
[240]1415              2(1X,E12.5),'Pa/m', &
1416             ' in x/y direction starting from dp_level_b =', F6.3, 'm', &
1417             A /)
[1]1418200 FORMAT (//' Run time and time step information:'/ &
1419             ' ----------------------------------'/)
1420201 FORMAT ( ' Timestep:          variable     maximum value: ',F6.3,' s', &
1421             '    CFL-factor: ',F4.2)
1422202 FORMAT ( ' Timestep:       dt = ',F6.3,' s'/)
1423203 FORMAT ( ' Start time:       ',F9.3,' s'/ &
1424             ' End time:         ',F9.3,' s')
1425204 FORMAT ( A,F9.3,' s')
1426205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
1427206 FORMAT (/' Time reached:     ',F9.3,' s'/ &
1428             ' CPU-time used:    ',F9.3,' s     per timestep:               ', &
1429               '  ',F9.3,' s'/                                                 &
1430             '                                   per second of simulated tim', &
1431               'e: ',F9.3,' s')
1432250 FORMAT (//' Computational grid and domain size:'/ &
1433              ' ----------------------------------'// &
1434              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1435              ' m    dz =    ',F7.3,' m'/ &
1436              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1437              ' m  z(u) = ',F10.3,' m'/)
1438252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
1439              ' factor: ',F5.3/ &
1440            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1441254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1442            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1443255 FORMAT (' Subdomains have equal size')
1444256 FORMAT (' Subdomains at the upper edges of the virtual processor grid ', &
1445              'have smaller sizes'/                                          &
1446            ' Size of smallest subdomain:    (  ',I4,',   ',I4,',   ',I4,')')
1447260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1448             ' degrees')
1449270 FORMAT (//' Topography informations:'/ &
1450              ' -----------------------'// &
1451              1X,'Topography: ',A)
1452271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1453              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1454                ' / ',I4)
[240]1455272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
1456              ' direction' / &
1457              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
1458              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
[138]1459280 FORMAT (//' Vegetation canopy (drag) model:'/ &
1460              ' ------------------------------'// &
1461              ' Canopy mode: ', A / &
1462              ' Canopy top: ',I4 / &
1463              ' Leaf drag coefficient: ',F6.2 /)
[153]1464281 FORMAT (/ ' Scalar_exchange_coefficient: ',F6.2 / &
1465              ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /)
1466282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s')
1467283 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
[138]1468              ' Height:              ',A,'  m'/ &
1469              ' Leaf area density:   ',A,'  m**2/m**3'/ &
1470              ' Gradient:            ',A,'  m**2/m**4'/ &
1471              ' Gridpoint:           ',A)
1472               
[1]1473300 FORMAT (//' Boundary conditions:'/ &
1474             ' -------------------'// &
1475             '                     p                    uv             ', &
1476             '                   pt'// &
1477             ' B. bound.: ',A/ &
1478             ' T. bound.: ',A)
[97]1479301 FORMAT (/'                     ',A// &
[1]1480             ' B. bound.: ',A/ &
1481             ' T. bound.: ',A)
[19]1482303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1483304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1484305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
1485               'computational u,v-level:'// &
[1]1486             '       zp = ',F6.2,' m   z0 = ',F6.4,' m   kappa = ',F4.2/ &
1487             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
[97]1488306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
[1]1489307 FORMAT ('       Heatflux has a random normal distribution')
1490308 FORMAT ('       Predefined surface temperature')
[97]1491309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
[1]1492310 FORMAT (//'    1D-Model:'// &
1493             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1494311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1495312 FORMAT ('       Predefined surface humidity')
1496313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1497314 FORMAT ('       Predefined scalar value at the surface')
[19]1498315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
[102]1499316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1500                    'atmosphere model')
[1]1501317 FORMAT (//' Lateral boundaries:'/ &
1502            '       left/right:  ',A/    &
1503            '       north/south: ',A)
1504318 FORMAT (/'       outflow damping layer width: ',I3,' gridpoints with km_', &
1505                    'max =',F5.1,' m**2/s')
[151]1506319 FORMAT ('       turbulence recycling at inflow switched on'/ &
1507            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1508            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
1509320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
[103]1510            '                                          v: ',F9.6,' m**2/s**2')
[151]1511325 FORMAT (//' List output:'/ &
[1]1512             ' -----------'//  &
1513            '    1D-Profiles:'/    &
1514            '       Output every             ',F8.2,' s')
[151]1515326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
[1]1516            '       Averaging input every    ',F8.2,' s')
1517330 FORMAT (//' Data output:'/ &
1518             ' -----------'/)
1519331 FORMAT (/'    1D-Profiles:')
1520332 FORMAT (/'       ',A)
1521333 FORMAT ('       Output every             ',F8.2,' s',/ &
1522            '       Time averaged over       ',F8.2,' s'/ &
1523            '       Averaging input every    ',F8.2,' s')
1524334 FORMAT (/'    2D-Arrays',A,':')
1525335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1526            '       Output every             ',F8.2,' s  ',A/ &
1527            '       Cross sections at ',A1,' = ',A/ &
1528            '       scalar-coordinates:   ',A,' m'/)
1529336 FORMAT (/'    3D-Arrays',A,':')
1530337 FORMAT (/'       Arrays: ',A/ &
1531            '       Output every             ',F8.2,' s  ',A/ &
1532            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1533338 FORMAT ('       Compressed data output'/ &
1534            '       Decimal precision: ',A/)
1535339 FORMAT ('       No output during initial ',F8.2,' s')
1536340 FORMAT (/'    Time series:')
1537341 FORMAT ('       Output every             ',F8.2,' s'/)
1538342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1539            '       Output every             ',F8.2,' s  ',A/ &
1540            '       Time averaged over       ',F8.2,' s'/ &
1541            '       Averaging input every    ',F8.2,' s'/ &
1542            '       Cross sections at ',A1,' = ',A/ &
1543            '       scalar-coordinates:   ',A,' m'/)
1544343 FORMAT (/'       Arrays: ',A/ &
1545            '       Output every             ',F8.2,' s  ',A/ &
1546            '       Time averaged over       ',F8.2,' s'/ &
1547            '       Averaging input every    ',F8.2,' s'/ &
1548            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1549345 FORMAT ('       Output format: ',A/)
1550#if defined( __dvrp_graphics )
1551360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
1552            '       Output every      ',F7.1,' s'/ &
1553            '       Output mode:      ',A/ &
1554            '       Host / User:      ',A,' / ',A/ &
1555            '       Directory:        ',A// &
1556            '       The sequence contains:')
1557361 FORMAT ('       Isosurface of ',A,'  Threshold value: ', E12.3)
1558362 FORMAT ('       Sectional plane ',A)
1559363 FORMAT ('       Particles')
[237]1560364 FORMAT (/'       Polygon reduction for topography: cluster_size = ', I1)
[1]1561#endif
1562#if defined( __spectra )
1563370 FORMAT ('    Spectra:')
1564371 FORMAT ('       Output every ',F7.1,' s'/)
1565372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
1566            '       Directions: ', 10(A5,',')/                         &
[189]1567            '       height levels  k = ', 20(I3,',')/                  &
1568            '                          ', 20(I3,',')/                  &
1569            '                          ', 20(I3,',')/                  &
1570            '                          ', 20(I3,',')/                  &
1571            '                          ', 19(I3,','),I3,'.'/           &
[1]1572            '       height levels selected for standard plot:'/        &
[189]1573            '                      k = ', 20(I3,',')/                  &
1574            '                          ', 20(I3,',')/                  &
1575            '                          ', 20(I3,',')/                  &
1576            '                          ', 20(I3,',')/                  &
1577            '                          ', 19(I3,','),I3,'.'/           &
[1]1578            '       Time averaged over ', F7.1, ' s,' /                &
1579            '       Profiles for the time averaging are taken every ', &
1580                    F6.1,' s')
1581#endif
1582400 FORMAT (//' Physical quantities:'/ &
1583              ' -------------------'/)
1584410 FORMAT ('    Angular velocity    :   omega = ',E9.3,' rad/s'/  &
1585            '    Geograph. latitude  :   phi   = ',F4.1,' degr'/   &
1586            '    Coriolis parameter  :   f     = ',F9.6,' 1/s'/    &
1587            '                            f*    = ',F9.6,' 1/s')
1588411 FORMAT (/'    Gravity             :   g     = ',F4.1,' m/s**2')
[97]1589412 FORMAT (/'    Reference density in buoyancy terms: ',F8.3,' kg/m**3')
1590413 FORMAT (/'    Reference temperature in buoyancy terms: ',F8.4,' K')
[57]1591415 FORMAT (/'    Cloud physics parameters:'/ &
[1]1592             '    ------------------------'/)
[57]1593416 FORMAT ('        Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
[1]1594            '        Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
1595            '        Density of air     :   rho_0 = ',F5.3,' kg/m**3'/  &
1596            '        Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
1597            '        Vapourization heat :   L_v   = ',E8.2,' J/kg')
1598420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1599            '       Height:        ',A,'  m'/ &
1600            '       Temperature:   ',A,'  K'/ &
1601            '       Gradient:      ',A,'  K/100m'/ &
1602            '       Gridpoint:     ',A)
1603421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1604            '       Height:      ',A,'  m'/ &
1605            '       Humidity:    ',A,'  kg/kg'/ &
1606            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1607            '       Gridpoint:   ',A)
1608422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1609            '       Height:                  ',A,'  m'/ &
1610            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1611            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1612            '       Gridpoint:               ',A)
1613423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1614            '       Height:      ',A,'  m'/ &
1615            '       ug:          ',A,'  m/s'/ &
1616            '       Gradient:    ',A,'  1/100s'/ &
1617            '       Gridpoint:   ',A)
1618424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1619            '       Height:      ',A,'  m'/ &
[97]1620            '       vg:          ',A,'  m/s'/ &
[1]1621            '       Gradient:    ',A,'  1/100s'/ &
1622            '       Gridpoint:   ',A)
[97]1623425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
1624            '       Height:     ',A,'  m'/ &
1625            '       Salinity:   ',A,'  psu'/ &
1626            '       Gradient:   ',A,'  psu/100m'/ &
1627            '       Gridpoint:  ',A)
[1]1628450 FORMAT (//' LES / Turbulence quantities:'/ &
1629              ' ---------------------------'/)
1630451 FORMAT ('   Diffusion coefficients are constant:'/ &
1631            '   Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1632452 FORMAT ('   Mixing length is limited to the Prandtl mixing lenth.')
1633453 FORMAT ('   Mixing length is limited to ',F4.2,' * z')
1634454 FORMAT ('   TKE is not allowed to fall below ',E9.2,' (m/s)**2')
[108]1635455 FORMAT ('   initial TKE is prescribed as ',E9.2,' (m/s)**2')
[1]1636470 FORMAT (//' Actions during the simulation:'/ &
1637              ' -----------------------------'/)
[94]1638471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1639            '    Disturbance amplitude           :     ',F4.2, ' m/s'/       &
1640            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1641            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
[1]1642472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1643                 ' to i/j =',I4)
1644473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1645                 1X,F5.3, ' m**2/s**2')
1646474 FORMAT ('    Random number generator used    : ',A/)
1647475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1648                 'respectively, if'/ &
1649            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1650                 ' 3D-simulation'/)
1651476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1652                 'respectively, if the'/ &
1653            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1654                 ' the 3D-simulation'/)
1655477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1656                 'respectively, if the'/ &
1657            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1658                 ' the 3D-simulation'/)
1659480 FORMAT ('    Particles:'/ &
1660            '    ---------'// &
1661            '       Particle advection is active (switched on at t = ', F7.1, &
1662                    ' s)'/ &
1663            '       Start of new particle generations every  ',F6.1,' s'/ &
1664            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
1665            '                            bottom:     ', A, ' top:         ', A/&
1666            '       Maximum particle age:                 ',F9.1,' s'/ &
[117]1667            '       Advection stopped at t = ',F9.1,' s'/ &
1668            '       Particles are sorted every ',F9.1,' s'/)
[1]1669481 FORMAT ('       Particles have random start positions'/)
1670482 FORMAT ('       Particles are advected only horizontally'/)
1671483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
1672484 FORMAT ('            Number of tails of the total domain: ',I10/ &
1673            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
1674            '            Maximum age of the end of the tail:  ',F8.2,' s')
1675485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
1676486 FORMAT ('       Particle statistics are written on file'/)
1677487 FORMAT ('       Number of particle groups: ',I2/)
1678488 FORMAT ('       SGS velocity components are used for particle advection'/ &
1679            '          minimum timestep for advection: ', F7.5/)
1680489 FORMAT ('       Number of particles simultaneously released at each ', &
1681                    'point: ', I5/)
1682490 FORMAT ('       Particle group ',I2,':'/ &
1683            '          Particle radius: ',E10.3, 'm')
1684491 FORMAT ('          Particle inertia is activated'/ &
1685            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
1686492 FORMAT ('          Particles are advected only passively (no inertia)'/)
1687493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
1688            '                                         y:',F8.1,' - ',F8.1,' m'/&
1689            '                                         z:',F8.1,' - ',F8.1,' m'/&
1690            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
1691                       ' m  dz = ',F8.1,' m'/)
1692494 FORMAT ('       Output of particle time series in NetCDF format every ', &
1693                    F8.2,' s'/)
1694495 FORMAT ('       Number of particles in total domain: ',I10/)
1695500 FORMAT (//' 1D-Model parameters:'/                           &
1696              ' -------------------'//                           &
1697            '    Simulation time:                   ',F8.1,' s'/ &
1698            '    Run-controll output every:         ',F8.1,' s'/ &
1699            '    Vertical profile output every:     ',F8.1,' s'/ &
1700            '    Mixing length calculation:         ',A/         &
1701            '    Dissipation calculation:           ',A/)
1702502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
1703
1704
1705 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.