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

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

preliminary uncomplete changes for ocean version

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