source: palm/tags/release-3.3/SOURCE/header.f90 @ 3890

Last change on this file since 3890 was 98, checked in by raasch, 17 years ago

updating comments and rc-file

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