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

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

preliminary version of modified boundary conditions at top

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