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

Last change on this file since 254 was 254, checked in by heinze, 15 years ago

Output of messages replaced by message handling routine.

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