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

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

vorlaeufige Standalone-Version fuer Linux-Cluster

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