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

Last change on this file since 156 was 153, checked in by steinfeld, 17 years ago

Update for the plant canopy model

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