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

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

Id keyword set as property for all *.f90 files

  • Property svn:keywords set to Id
File size: 53.8 KB
Line 
1 SUBROUTINE header
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6!
7!
8! Former revisions:
9! -----------------
10! $Id: header.f90 4 2007-02-13 11:33:16Z 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) = pt_top'
348    ELSE
349       roben  = TRIM( roben  ) // ' pt(nzt) = pt(nzt-1) + dpt/dz'
350    ENDIF
351
352    WRITE ( io, 300 )  runten, roben
353
354    IF ( .NOT. constant_diffusion )  THEN
355       IF ( ibc_e_b == 1 )  THEN
356          runten = 'e(0)     = e(1)'
357       ELSE
358          runten = 'e(0)     = e(1) = (u*/0.1)**2'
359       ENDIF
360       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
361
362       WRITE ( io, 301 )  runten, roben       
363
364    ENDIF
365
366    IF ( moisture  .OR.  passive_scalar )  THEN
367       IF ( moisture )  THEN
368          IF ( ibc_q_b == 0 )  THEN
369             runten = 'q(0)     = q_surface'
370          ELSE
371             runten = 'q(0)     = q(1)'
372          ENDIF
373          IF ( ibc_q_t == 0 )  THEN
374             roben =  'q(nzt)   = q_top'
375          ELSE
376             roben =  'q(nzt)   = q(nzt-1) + dq/dz'
377          ENDIF
378       ELSE
379          IF ( ibc_q_b == 0 )  THEN
380             runten = 's(0)     = s_surface'
381          ELSE
382             runten = 's(0)     = s(1)'
383          ENDIF
384          IF ( ibc_q_t == 0 )  THEN
385             roben =  's(nzt)   = s_top'
386          ELSE
387             roben =  's(nzt)   = s(nzt-1) + ds/dz'
388          ENDIF
389       ENDIF
390
391       WRITE ( io, 302 ) runten, roben
392
393    ENDIF
394
395    IF ( use_surface_fluxes )  THEN
396       WRITE ( io, 303 )
397       IF ( constant_heatflux )  THEN
398          WRITE ( io, 306 )  surface_heatflux
399          IF ( random_heatflux )  WRITE ( io, 307 )
400       ENDIF
401       IF ( moisture  .AND.  constant_waterflux )  THEN
402          WRITE ( io, 311 ) surface_waterflux
403       ENDIF
404       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
405          WRITE ( io, 313 ) surface_waterflux
406       ENDIF
407    ENDIF
408
409    IF ( prandtl_layer )  THEN
410       WRITE ( io, 305 )  zu(1), roughness_length, kappa, rif_min, rif_max
411       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
412       IF ( moisture  .AND.  .NOT. constant_waterflux )  THEN
413          WRITE ( io, 312 )
414       ENDIF
415       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
416          WRITE ( io, 314 )
417       ENDIF
418    ELSE
419       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
420          WRITE ( io, 310 )  rif_min, rif_max
421       ENDIF
422    ENDIF
423
424    WRITE ( io, 317 )  bc_lr, bc_ns
425    IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
426       WRITE ( io, 318 )  outflow_damping_width, km_damp_max
427    ENDIF
428
429!
430!-- Listing of 1D-profiles
431    WRITE ( io, 320 )  dt_dopr_listing
432    IF ( averaging_interval_pr /= 0.0 )  THEN
433       WRITE ( io, 321 )  averaging_interval_pr, dt_averaging_input_pr
434    ENDIF
435
436!
437!-- DATA output
438    WRITE ( io, 330 )
439    IF ( averaging_interval_pr /= 0.0 )  THEN
440       WRITE ( io, 321 )  averaging_interval_pr, dt_averaging_input_pr
441    ENDIF
442
443!
444!-- 1D-profiles
445    dopr_chr = 'Profile:'
446    IF ( dopr_n /= 0 )  THEN
447       WRITE ( io, 331 )
448
449       output_format = ''
450       IF ( netcdf_output )  THEN
451          IF ( netcdf_64bit )  THEN
452             output_format = 'netcdf (64 bit offset)'
453          ELSE
454             output_format = 'netcdf'
455          ENDIF
456       ENDIF
457       IF ( profil_output )  THEN
458          IF ( netcdf_output )  THEN
459             output_format = TRIM( output_format ) // ' and profil'
460          ELSE
461             output_format = 'profil'
462          ENDIF
463       ENDIF
464       WRITE ( io, 345 )  output_format
465
466       DO  i = 1, dopr_n
467          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
468          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
469             WRITE ( io, 332 )  dopr_chr
470             dopr_chr = '       :'
471          ENDIF
472       ENDDO
473
474       IF ( dopr_chr /= '' )  THEN
475          WRITE ( io, 332 )  dopr_chr
476       ENDIF
477       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
478       IF ( skip_time_dopr /= 0.0 )  WRITE ( io, 339 )  skip_time_dopr
479    ENDIF
480
481!
482!-- 2D-arrays
483    DO  av = 0, 1
484
485       i = 1
486       do2d_xy = ''
487       do2d_xz = ''
488       do2d_yz = ''
489       DO  WHILE ( do2d(av,i) /= ' ' )
490
491          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
492          do2d_mode = do2d(av,i)(l-1:l)
493
494          SELECT CASE ( do2d_mode )
495             CASE ( 'xy' )
496                ll = LEN_TRIM( do2d_xy )
497                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
498             CASE ( 'xz' )
499                ll = LEN_TRIM( do2d_xz )
500                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
501             CASE ( 'yz' )
502                ll = LEN_TRIM( do2d_yz )
503                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
504          END SELECT
505
506          i = i + 1
507
508       ENDDO
509
510       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
511              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
512              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) )  .AND. &
513            ( netcdf_output  .OR.  iso2d_output ) )  THEN
514
515          IF (  av == 0 )  THEN
516             WRITE ( io, 334 )  ''
517          ELSE
518             WRITE ( io, 334 )  '(time-averaged)'
519          ENDIF
520
521          IF ( do2d_at_begin )  THEN
522             begin_chr = 'and at the start'
523          ELSE
524             begin_chr = ''
525          ENDIF
526
527          output_format = ''
528          IF ( netcdf_output )  THEN
529             IF ( netcdf_64bit )  THEN
530                output_format = 'netcdf (64 bit offset)'
531             ELSE
532                output_format = 'netcdf'
533             ENDIF
534          ENDIF
535          IF ( iso2d_output )  THEN
536             IF ( netcdf_output )  THEN
537                output_format = TRIM( output_format ) // ' and iso2d'
538             ELSE
539                output_format = 'iso2d'
540             ENDIF
541          ENDIF
542          WRITE ( io, 345 )  output_format
543
544          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
545             i = 1
546             slices = '/'
547             coordinates = '/'
548!
549!--          Building strings with index and coordinate informations of the
550!--          slices
551             DO  WHILE ( section(i,1) /= -9999 )
552
553                WRITE (section_chr,'(I5)')  section(i,1)
554                section_chr = ADJUSTL( section_chr )
555                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
556
557                WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
558                coor_chr = ADJUSTL( coor_chr )
559                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
560
561                i = i + 1
562             ENDDO
563             IF ( av == 0 )  THEN
564                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
565                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
566                                   TRIM( coordinates )
567                IF ( skip_time_do2d_xy /= 0.0 )  THEN
568                   WRITE ( io, 339 )  skip_time_do2d_xy
569                ENDIF
570             ELSE
571                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
572                                   TRIM( begin_chr ), averaging_interval, &
573                                   dt_averaging_input, 'k', TRIM( slices ), &
574                                   TRIM( coordinates )
575                IF ( skip_time_data_output_av /= 0.0 )  THEN
576                   WRITE ( io, 339 )  skip_time_data_output_av
577                ENDIF
578             ENDIF
579
580          ENDIF
581
582          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
583             i = 1
584             slices = '/'
585             coordinates = '/'
586!
587!--          Building strings with index and coordinate informations of the
588!--          slices
589             DO  WHILE ( section(i,2) /= -9999 )
590
591                WRITE (section_chr,'(I5)')  section(i,2)
592                section_chr = ADJUSTL( section_chr )
593                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
594
595                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
596                coor_chr = ADJUSTL( coor_chr )
597                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
598
599                i = i + 1
600             ENDDO
601             IF ( av == 0 )  THEN
602                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
603                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
604                                   TRIM( coordinates )
605                IF ( skip_time_do2d_xz /= 0.0 )  THEN
606                   WRITE ( io, 339 )  skip_time_do2d_xz
607                ENDIF
608             ELSE
609                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
610                                   TRIM( begin_chr ), averaging_interval, &
611                                   dt_averaging_input, 'j', TRIM( slices ), &
612                                   TRIM( coordinates )
613                IF ( skip_time_data_output_av /= 0.0 )  THEN
614                   WRITE ( io, 339 )  skip_time_data_output_av
615                ENDIF
616             ENDIF
617          ENDIF
618
619          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
620             i = 1
621             slices = '/'
622             coordinates = '/'
623!
624!--          Building strings with index and coordinate informations of the
625!--          slices
626             DO  WHILE ( section(i,3) /= -9999 )
627
628                WRITE (section_chr,'(I5)')  section(i,3)
629                section_chr = ADJUSTL( section_chr )
630                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
631
632                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
633                coor_chr = ADJUSTL( coor_chr )
634                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
635
636                i = i + 1
637             ENDDO
638             IF ( av == 0 )  THEN
639                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
640                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
641                                   TRIM( coordinates )
642                IF ( skip_time_do2d_yz /= 0.0 )  THEN
643                   WRITE ( io, 339 )  skip_time_do2d_yz
644                ENDIF
645             ELSE
646                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
647                                   TRIM( begin_chr ), averaging_interval, &
648                                   dt_averaging_input, 'i', TRIM( slices ), &
649                                   TRIM( coordinates )
650                IF ( skip_time_data_output_av /= 0.0 )  THEN
651                   WRITE ( io, 339 )  skip_time_data_output_av
652                ENDIF
653             ENDIF
654          ENDIF
655
656       ENDIF
657
658    ENDDO
659
660!
661!-- 3d-arrays
662    DO  av = 0, 1
663
664       i = 1
665       do3d_chr = ''
666       DO  WHILE ( do3d(av,i) /= ' ' )
667
668          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
669          i = i + 1
670
671       ENDDO
672
673       IF ( do3d_chr /= '' )  THEN
674          IF ( av == 0 )  THEN
675             WRITE ( io, 336 )  ''
676          ELSE
677             WRITE ( io, 336 )  '(time-averaged)'
678          ENDIF
679
680          output_format = ''
681          IF ( netcdf_output )  THEN
682             IF ( netcdf_64bit )  THEN
683                output_format = 'netcdf (64 bit offset)'
684             ELSE
685                output_format = 'netcdf'
686             ENDIF
687          ENDIF
688          IF ( avs_output )  THEN
689             IF ( netcdf_output )  THEN
690                output_format = TRIM( output_format ) // ' and avs'
691             ELSE
692                output_format = 'avs'
693             ENDIF
694          ENDIF
695          WRITE ( io, 345 )  output_format
696
697          IF ( do3d_at_begin )  THEN
698             begin_chr = 'and at the start'
699          ELSE
700             begin_chr = ''
701          ENDIF
702          IF ( av == 0 )  THEN
703             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
704                                zu(nz_do3d), nz_do3d
705          ELSE
706             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
707                                TRIM( begin_chr ), averaging_interval, &
708                                dt_averaging_input, zu(nz_do3d), nz_do3d
709          ENDIF
710
711          IF ( do3d_compress )  THEN
712             do3d_chr = ''
713             i = 1
714             DO WHILE ( do3d(av,i) /= ' ' )
715
716                SELECT CASE ( do3d(av,i) )
717                   CASE ( 'u' )
718                      j = 1
719                   CASE ( 'v' )
720                      j = 2
721                   CASE ( 'w' )
722                      j = 3
723                   CASE ( 'p' )
724                      j = 4
725                   CASE ( 'pt' )
726                      j = 5
727                END SELECT
728                WRITE ( prec, '(I1)' )  plot_3d_precision(j)%precision
729                do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // &
730                           ':' // prec // ','
731                i = i + 1
732
733             ENDDO
734             WRITE ( io, 338 )  do3d_chr
735
736          ENDIF
737
738          IF ( av == 0 )  THEN
739             IF ( skip_time_do3d /= 0.0 )  THEN
740                WRITE ( io, 339 )  skip_time_do3d
741             ENDIF
742          ELSE
743             IF ( skip_time_data_output_av /= 0.0 )  THEN
744                WRITE ( io, 339 )  skip_time_data_output_av
745             ENDIF
746          ENDIF
747
748       ENDIF
749
750    ENDDO
751
752!
753!-- Timeseries
754    IF ( dt_dots /= 9999999.9 )  THEN
755       WRITE ( io, 340 )
756
757       output_format = ''
758       IF ( netcdf_output )  THEN
759          IF ( netcdf_64bit )  THEN
760             output_format = 'netcdf (64 bit offset)'
761          ELSE
762             output_format = 'netcdf'
763          ENDIF
764       ENDIF
765       IF ( profil_output )  THEN
766          IF ( netcdf_output )  THEN
767             output_format = TRIM( output_format ) // ' and profil'
768          ELSE
769             output_format = 'profil'
770          ENDIF
771       ENDIF
772       WRITE ( io, 345 )  output_format
773       WRITE ( io, 341 )  dt_dots
774    ENDIF
775
776#if defined( __dvrp_graphics )
777!
778!-- Dvrp-output
779    IF ( dt_dvrp /= 9999999.9 )  THEN
780       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
781                          TRIM( dvrp_username ), TRIM( dvrp_directory )
782       i = 1
783       l = 0
784       DO WHILE ( mode_dvrp(i) /= ' ' )
785          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
786             READ ( mode_dvrp(i), '(10X,I1)' )  j
787             l = l + 1
788             IF ( do3d(0,j) /= ' ' )  THEN
789                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l)
790             ENDIF
791          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
792             READ ( mode_dvrp(i), '(6X,I1)' )  j
793             IF ( do2d(0,j) /= ' ' )  WRITE ( io, 362 )  TRIM( do2d(0,j) )
794          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
795             WRITE ( io, 363 )
796          ENDIF
797          i = i + 1
798       ENDDO
799    ENDIF
800#endif
801
802#if defined( __spectra )
803!
804!-- Spectra output
805    IF ( dt_dosp /= 9999999.9 ) THEN
806       WRITE ( io, 370 )
807
808       output_format = ''
809       IF ( netcdf_output )  THEN
810          IF ( netcdf_64bit )  THEN
811             output_format = 'netcdf (64 bit offset)'
812          ELSE
813             output_format = 'netcdf'
814          ENDIF
815       ENDIF
816       IF ( profil_output )  THEN
817          IF ( netcdf_output )  THEN
818             output_format = TRIM( output_format ) // ' and profil'
819          ELSE
820             output_format = 'profil'
821          ENDIF
822       ENDIF
823       WRITE ( io, 345 )  output_format
824       WRITE ( io, 371 )  dt_dosp
825       IF ( skip_time_dosp /= 0.0 )  WRITE ( io, 339 )  skip_time_dosp
826       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
827                          ( spectra_direction(i), i = 1,10 ),  &
828                          ( comp_spectra_level(i), i = 1,10 ), &
829                          ( plot_spectra_level(i), i = 1,10 ), &
830                          averaging_interval_sp, dt_averaging_input_pr
831    ENDIF
832#endif
833
834    WRITE ( io, 99 )
835
836!
837!-- Physical quantities
838    WRITE ( io, 400 )
839
840!
841!-- Geostrophic parameters
842    WRITE ( io, 410 )  omega, phi, f, fs
843
844!
845!-- Other quantities
846    WRITE ( io, 411 )  g
847
848!
849!-- Cloud physics parameters
850    IF ( cloud_physics ) THEN
851       WRITE ( io, 412 )
852       WRITE ( io, 413 ) surface_pressure, r_d, rho_surface, cp, l_v
853    ENDIF
854
855!-- Profile of the geostrophic wind (component ug)
856!-- Building output strings
857    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
858    gradients = '------'
859    slices = '     0'
860    coordinates = '   0.0'
861    i = 1
862    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
863     
864       WRITE (coor_chr,'(F6.2,4X)')  ug(ug_vertical_gradient_level_ind(i))
865       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
866
867       WRITE (coor_chr,'(F6.2,4X)')  ug_vertical_gradient(i)
868       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
869
870       WRITE (coor_chr,'(I6,4X)')  ug_vertical_gradient_level_ind(i)
871       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
872
873       WRITE (coor_chr,'(F6.1,4X)')  ug_vertical_gradient_level(i)
874       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
875
876       i = i + 1
877    ENDDO
878
879    WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
880                       TRIM( gradients ), TRIM( slices )
881
882!-- Profile of the geostrophic wind (component vg)
883!-- Building output strings
884    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
885    gradients = '------'
886    slices = '     0'
887    coordinates = '   0.0'
888    i = 1
889    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
890
891       WRITE (coor_chr,'(F6.2,4X)')  vg(vg_vertical_gradient_level_ind(i))
892       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
893
894       WRITE (coor_chr,'(F6.2,4X)')  vg_vertical_gradient(i)
895       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
896
897       WRITE (coor_chr,'(I6,4X)')  vg_vertical_gradient_level_ind(i)
898       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
899
900       WRITE (coor_chr,'(F6.1,4X)')  vg_vertical_gradient_level(i)
901       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
902
903       i = i + 1 
904    ENDDO
905
906    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
907                       TRIM( gradients ), TRIM( slices )
908
909!
910!-- Initial temperature profile
911!-- Building output strings, starting with surface temperature
912    WRITE ( temperatures, '(F6.2)' )  pt_surface
913    gradients = '------'
914    slices = '     0'
915    coordinates = '   0.0'
916    i = 1
917    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
918
919       WRITE (coor_chr,'(F6.2,4X)')  pt_init(pt_vertical_gradient_level_ind(i))
920       temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
921
922       WRITE (coor_chr,'(F6.2,4X)')  pt_vertical_gradient(i)
923       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
924
925       WRITE (coor_chr,'(I6,4X)')  pt_vertical_gradient_level_ind(i)
926       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
927
928       WRITE (coor_chr,'(F6.1,4X)')  pt_vertical_gradient_level(i)
929       coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
930
931       i = i + 1
932    ENDDO
933
934    WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
935                       TRIM( gradients ), TRIM( slices )
936
937!
938!-- Initial humidity profile
939!-- Building output strings, starting with surface humidity
940    IF ( moisture  .OR.  passive_scalar )  THEN
941       WRITE ( temperatures, '(E8.1)' )  q_surface
942       gradients = '--------'
943       slices = '       0'
944       coordinates = '     0.0'
945       i = 1
946       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
947         
948          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
949          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
950
951          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
952          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
953         
954          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
955          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
956         
957          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
958          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
959
960          i = i + 1
961       ENDDO
962
963       IF ( moisture )  THEN
964          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
965                             TRIM( gradients ), TRIM( slices )
966       ELSE
967          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
968                             TRIM( gradients ), TRIM( slices )
969       ENDIF
970    ENDIF
971
972!
973!-- LES / turbulence parameters
974    WRITE ( io, 450 )
975
976!--
977! ... LES-constants used must still be added here
978!--
979    IF ( constant_diffusion )  THEN
980       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
981                          prandtl_number
982    ENDIF
983    IF ( .NOT. constant_diffusion)  THEN
984       IF ( e_min > 0.0 )  WRITE ( io, 454 )  e_min
985       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
986       IF ( adjust_mixing_length  .AND.  prandtl_layer )  WRITE ( io, 452 )
987    ENDIF
988
989!
990!-- Special actions during the run
991    WRITE ( io, 470 )
992    IF ( create_disturbances )  THEN
993       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
994                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
995                          zu(disturbance_level_ind_t), disturbance_level_ind_t
996       IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
997          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
998       ELSE
999          WRITE ( io, 473 )  disturbance_energy_limit
1000       ENDIF
1001       WRITE ( io, 474 )  TRIM( random_generator )
1002    ENDIF
1003    IF ( pt_surface_initial_change /= 0.0 )  THEN
1004       WRITE ( io, 475 )  pt_surface_initial_change
1005    ENDIF
1006    IF ( moisture  .AND.  q_surface_initial_change /= 0.0 )  THEN
1007       WRITE ( io, 476 )  q_surface_initial_change       
1008    ENDIF
1009    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0 )  THEN
1010       WRITE ( io, 477 )  q_surface_initial_change       
1011    ENDIF
1012
1013#if defined( __particles )
1014!
1015!-- Particle attributes
1016    WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1017                       bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
1018                       end_time_prel
1019    IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1020    IF ( random_start_position )  WRITE ( io, 481 )
1021    IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1022    WRITE ( io, 495 )  total_number_of_particles
1023    IF ( .NOT. vertical_particle_advection )  WRITE ( io, 482 )
1024    IF ( maximum_number_of_tailpoints /= 0 )  THEN
1025       WRITE ( io, 483 )  maximum_number_of_tailpoints
1026       IF ( minimum_tailpoint_distance /= 0 )  THEN
1027          WRITE ( io, 484 )  total_number_of_tails, minimum_tailpoint_distance,&
1028                             maximum_tailpoint_age
1029       ENDIF
1030    ENDIF
1031    IF ( dt_write_particle_data /= 9999999.9 )  THEN
1032       WRITE ( io, 485 )  dt_write_particle_data
1033       output_format = ''
1034       IF ( netcdf_output )  THEN
1035          IF ( netcdf_64bit )  THEN
1036             output_format = 'netcdf (64 bit offset) and binary'
1037          ELSE
1038             output_format = 'netcdf and binary'
1039          ENDIF
1040       ELSE
1041          output_format = 'binary'
1042       ENDIF
1043       WRITE ( io, 345 )  output_format
1044    ENDIF
1045    IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1046    IF ( write_particle_statistics )  WRITE ( io, 486 )
1047
1048    WRITE ( io, 487 )  number_of_particle_groups
1049
1050    DO  i = 1, number_of_particle_groups
1051       IF ( i == 1  .AND.  density_ratio(i) == 9999999.9 )  THEN
1052          WRITE ( io, 490 )  i, 0.0
1053          WRITE ( io, 492 )
1054       ELSE
1055          WRITE ( io, 490 )  i, radius(i)
1056          IF ( density_ratio(i) /= 0.0 )  THEN
1057             WRITE ( io, 491 )  density_ratio(i)
1058          ELSE
1059             WRITE ( io, 492 )
1060          ENDIF
1061       ENDIF
1062       WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1063                          pdx(i), pdy(i), pdz(i)
1064    ENDDO
1065
1066#endif
1067
1068!
1069!-- Parameters of 1D-model
1070    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1071       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1072                          mixing_length_1d, dissipation_1d
1073       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1074          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1075       ENDIF
1076    ENDIF
1077
1078!
1079!-- User-defined informations
1080    CALL user_header( io )
1081
1082    WRITE ( io, 99 )
1083
1084#if defined( __ibm )
1085!
1086!-- Write buffer contents to disc immediately
1087    CALL FLUSH_( io )
1088#elif defined( __lcmuk )  ||  defined( __nec )
1089    CALL FLUSH( io )
1090#endif
1091
1092!
1093!-- Here the FORMATs start
1094
1095 99 FORMAT (1X,78('-'))
1096100 FORMAT (/10X,'****************',11X,28('-')/                &
1097            10X,'*  ',A12,'*',11X,A/                            &
1098            10X,'****************',11X,28('-')//                &
1099            ' Date:            ',A8,11X,'Run:       ',A20/      &
1100            ' Time:            ',A8,11X,'Run-No.:   ',I2.2/     &
1101            ' Run on host:   ',A10)
1102#if defined( __parallel )
1103101 FORMAT (' Number of PEs:',7X,I4,11X,'Processor grid (x,y): (',I3,',',I3, &
1104              ')',1X,A)
1105102 FORMAT (' Number of PEs:',7X,I4,11X,'Tasks:',I4,'   threads per task:',I4/ &
1106              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
1107103 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1108104 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
1109            37X,'because the job is running on an SMP-cluster')
1110105 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
1111#endif
1112110 FORMAT (/' Numerical Schemes:'/ &
1113             ' -----------------'/)
1114111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1115112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1116            '     Iterations (initial/other): ',I3,'/',I3,'  omega = ',F5.3)
1117113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1118                  ' or Upstream')
1119114 FORMAT (' --> Momentum advection via Upstream-Spline-Scheme')
1120115 FORMAT ('     Tendencies are smoothed via Long-Filter with factor ',F5.3) 
1121116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1122                  ' or Upstream')
1123117 FORMAT (' --> Scalar advection via Upstream-Spline-Scheme')
1124118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1125119 FORMAT (' --> Galilei-Transform applied to horizontal advection', &
1126            '     Translation velocity = ',A/ &
1127            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1128120 FORMAT (' --> Time differencing scheme: leapfrog only (no euler in case', &
1129                  ' of timestep changes)')
1130121 FORMAT (' --> Time differencing scheme: leapfrog + euler in case of', &
1131                  ' timestep changes')
1132122 FORMAT (' --> Time differencing scheme: ',A)
1133123 FORMAT (' --> Rayleigh-Damping active, starts above z = ',F8.2,' m'/ &
1134            '     maximum damping coefficient: ',F5.3, ' 1/s')
1135124 FORMAT ('     Spline-overshoots are being suppressed')
1136125 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1137                  ' of'/                                                       &
1138            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1139126 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1140                  ' of'/                                                       &
1141            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1142127 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1143            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1144128 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1145            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1146129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1147130 FORMAT (' --> Additional prognostic equation for the total water content')
1148131 FORMAT (' --> Parameterization of condensation processes via (0%-or100%)')
1149132 FORMAT (' --> Parameterization of long-wave radiation processes via'/ &
1150            '     effective emissivity scheme')
1151133 FORMAT (' --> Precipitation parameterization via Kessler-Scheme')
1152134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1153135 FORMAT (' --> Solve perturbation pressure via multigrid method (', &
1154                  A,'-cycle)'/ &
1155            '     number of grid levels:                   ',I2/ &
1156            '     Gauss-Seidel red/black iterations:       ',I2)
1157136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1158                  I3,')')
1159137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1160            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1161                  I3,')'/ &
1162            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1163                  I3,')')
1164138 FORMAT ('     Using hybrid version for 1d-domain-decomposition')
1165139 FORMAT (' --> Prognostic equations are solved in one single loop (fast', &
1166                  ' method)')
1167140 FORMAT ('     maximum residual allowed:                ',E10.3)
1168141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1169142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1170                  'step')
1171150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
1172                  'conserved')
1173200 FORMAT (//' Run time and time step information:'/ &
1174             ' ----------------------------------'/)
1175201 FORMAT ( ' Timestep:          variable     maximum value: ',F6.3,' s', &
1176             '    CFL-factor: ',F4.2)
1177202 FORMAT ( ' Timestep:       dt = ',F6.3,' s'/)
1178203 FORMAT ( ' Start time:       ',F9.3,' s'/ &
1179             ' End time:         ',F9.3,' s')
1180204 FORMAT ( A,F9.3,' s')
1181205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
1182206 FORMAT (/' Time reached:     ',F9.3,' s'/ &
1183             ' CPU-time used:    ',F9.3,' s     per timestep:               ', &
1184               '  ',F9.3,' s'/                                                 &
1185             '                                   per second of simulated tim', &
1186               'e: ',F9.3,' s')
1187250 FORMAT (//' Computational grid and domain size:'/ &
1188              ' ----------------------------------'// &
1189              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1190              ' m    dz =    ',F7.3,' m'/ &
1191              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1192              ' m  z(u) = ',F10.3,' m'/)
1193252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
1194              ' factor: ',F5.3/ &
1195            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1196254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1197            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1198255 FORMAT (' Subdomains have equal size')
1199256 FORMAT (' Subdomains at the upper edges of the virtual processor grid ', &
1200              'have smaller sizes'/                                          &
1201            ' Size of smallest subdomain:    (  ',I4,',   ',I4,',   ',I4,')')
1202260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1203             ' degrees')
1204270 FORMAT (//' Topography informations:'/ &
1205              ' -----------------------'// &
1206              1X,'Topography: ',A)
1207271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1208              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1209                ' / ',I4)
1210300 FORMAT (//' Boundary conditions:'/ &
1211             ' -------------------'// &
1212             '                     p                    uv             ', &
1213             '                   pt'// &
1214             ' B. bound.: ',A/ &
1215             ' T. bound.: ',A)
1216301 FORMAT (/'                     e'// &
1217             ' B. bound.: ',A/ &
1218             ' T. bound.: ',A)
1219302 FORMAT (/'                     q'// &
1220             ' B. bound.: ',A/ &
1221             ' T. bound.: ',A)
1222303 FORMAT (/' Surface fluxes are used in diffusion terms at k=1')
1223305 FORMAT (//'    Prandtl-Layer between surface and first computational ', &
1224               'u,v-level:'// &
1225             '       zp = ',F6.2,' m   z0 = ',F6.4,' m   kappa = ',F4.2/ &
1226             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1227306 FORMAT ('       Predefined constant heatflux:   ',F6.3,' K m/s')
1228307 FORMAT ('       Heatflux has a random normal distribution')
1229308 FORMAT ('       Predefined surface temperature')
1230310 FORMAT (//'    1D-Model:'// &
1231             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1232311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1233312 FORMAT ('       Predefined surface humidity')
1234313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1235314 FORMAT ('       Predefined scalar value at the surface')
1236317 FORMAT (//' Lateral boundaries:'/ &
1237            '       left/right:  ',A/    &
1238            '       north/south: ',A)
1239318 FORMAT (/'       outflow damping layer width: ',I3,' gridpoints with km_', &
1240                    'max =',F5.1,' m**2/s')
1241320 FORMAT (//' List output:'/ &
1242             ' -----------'//  &
1243            '    1D-Profiles:'/    &
1244            '       Output every             ',F8.2,' s')
1245321 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
1246            '       Averaging input every    ',F8.2,' s')
1247330 FORMAT (//' Data output:'/ &
1248             ' -----------'/)
1249331 FORMAT (/'    1D-Profiles:')
1250332 FORMAT (/'       ',A)
1251333 FORMAT ('       Output every             ',F8.2,' s',/ &
1252            '       Time averaged over       ',F8.2,' s'/ &
1253            '       Averaging input every    ',F8.2,' s')
1254334 FORMAT (/'    2D-Arrays',A,':')
1255335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1256            '       Output every             ',F8.2,' s  ',A/ &
1257            '       Cross sections at ',A1,' = ',A/ &
1258            '       scalar-coordinates:   ',A,' m'/)
1259336 FORMAT (/'    3D-Arrays',A,':')
1260337 FORMAT (/'       Arrays: ',A/ &
1261            '       Output every             ',F8.2,' s  ',A/ &
1262            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1263338 FORMAT ('       Compressed data output'/ &
1264            '       Decimal precision: ',A/)
1265339 FORMAT ('       No output during initial ',F8.2,' s')
1266340 FORMAT (/'    Time series:')
1267341 FORMAT ('       Output every             ',F8.2,' s'/)
1268342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1269            '       Output every             ',F8.2,' s  ',A/ &
1270            '       Time averaged over       ',F8.2,' s'/ &
1271            '       Averaging input every    ',F8.2,' s'/ &
1272            '       Cross sections at ',A1,' = ',A/ &
1273            '       scalar-coordinates:   ',A,' m'/)
1274343 FORMAT (/'       Arrays: ',A/ &
1275            '       Output every             ',F8.2,' s  ',A/ &
1276            '       Time averaged over       ',F8.2,' s'/ &
1277            '       Averaging input every    ',F8.2,' s'/ &
1278            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1279345 FORMAT ('       Output format: ',A/)
1280#if defined( __dvrp_graphics )
1281360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
1282            '       Output every      ',F7.1,' s'/ &
1283            '       Output mode:      ',A/ &
1284            '       Host / User:      ',A,' / ',A/ &
1285            '       Directory:        ',A// &
1286            '       The sequence contains:')
1287361 FORMAT ('       Isosurface of ',A,'  Threshold value: ', E12.3)
1288362 FORMAT ('       Sectional plane ',A)
1289363 FORMAT ('       Particles')
1290#endif
1291#if defined( __spectra )
1292370 FORMAT ('    Spectra:')
1293371 FORMAT ('       Output every ',F7.1,' s'/)
1294372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
1295            '       Directions: ', 10(A5,',')/                         &
1296            '       height levels  k = ', 9(I3,','),I3,'.'/            &
1297            '       height levels selected for standard plot:'/        &
1298            '                      k = ', 9(I3,','),I3,'.'/            &
1299            '       Time averaged over ', F7.1, ' s,' /                &
1300            '       Profiles for the time averaging are taken every ', &
1301                    F6.1,' s')
1302#endif
1303400 FORMAT (//' Physical quantities:'/ &
1304              ' -------------------'/)
1305410 FORMAT ('    Angular velocity    :   omega = ',E9.3,' rad/s'/  &
1306            '    Geograph. latitude  :   phi   = ',F4.1,' degr'/   &
1307            '    Coriolis parameter  :   f     = ',F9.6,' 1/s'/    &
1308            '                            f*    = ',F9.6,' 1/s')
1309411 FORMAT (/'    Gravity             :   g     = ',F4.1,' m/s**2')
1310412 FORMAT (/'    Cloud physics parameters:'/ &
1311             '    ------------------------'/)
1312413 FORMAT ('        Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
1313            '        Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
1314            '        Density of air     :   rho_0 = ',F5.3,' kg/m**3'/  &
1315            '        Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
1316            '        Vapourization heat :   L_v   = ',E8.2,' J/kg')
1317420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1318            '       Height:        ',A,'  m'/ &
1319            '       Temperature:   ',A,'  K'/ &
1320            '       Gradient:      ',A,'  K/100m'/ &
1321            '       Gridpoint:     ',A)
1322421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1323            '       Height:      ',A,'  m'/ &
1324            '       Humidity:    ',A,'  kg/kg'/ &
1325            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1326            '       Gridpoint:   ',A)
1327422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1328            '       Height:                  ',A,'  m'/ &
1329            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1330            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1331            '       Gridpoint:               ',A)
1332423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1333            '       Height:      ',A,'  m'/ &
1334            '       ug:          ',A,'  m/s'/ &
1335            '       Gradient:    ',A,'  1/100s'/ &
1336            '       Gridpoint:   ',A)
1337424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1338            '       Height:      ',A,'  m'/ &
1339            '       vg:          ',A,'  m/S'/ &
1340            '       Gradient:    ',A,'  1/100s'/ &
1341            '       Gridpoint:   ',A)
1342450 FORMAT (//' LES / Turbulence quantities:'/ &
1343              ' ---------------------------'/)
1344451 FORMAT ('   Diffusion coefficients are constant:'/ &
1345            '   Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1346452 FORMAT ('   Mixing length is limited to the Prandtl mixing lenth.')
1347453 FORMAT ('   Mixing length is limited to ',F4.2,' * z')
1348454 FORMAT ('   TKE is not allowed to fall below ',E9.2,' (m/s)**2')
1349470 FORMAT (//' Actions during the simulation:'/ &
1350              ' -----------------------------'/)
1351471 FORMAT ('    Disturbance impulse (u,v) every :  ',F6.2,' s'/             &
1352            '    Disturbance amplitude           :    ',F4.2, ' m/s'/        &
1353            '    Lower disturbance level         : ',F7.2,' m (GP ',I4,')'/  &
1354            '    Upper disturbance level         : ',F7.2,' m (GP ',I4,')')
1355472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1356                 ' to i/j =',I4)
1357473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1358                 1X,F5.3, ' m**2/s**2')
1359474 FORMAT ('    Random number generator used    : ',A/)
1360475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1361                 'respectively, if'/ &
1362            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1363                 ' 3D-simulation'/)
1364476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1365                 'respectively, if the'/ &
1366            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1367                 ' the 3D-simulation'/)
1368477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1369                 'respectively, if the'/ &
1370            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1371                 ' the 3D-simulation'/)
1372#if defined( __particles )
1373480 FORMAT ('    Particles:'/ &
1374            '    ---------'// &
1375            '       Particle advection is active (switched on at t = ', F7.1, &
1376                    ' s)'/ &
1377            '       Start of new particle generations every  ',F6.1,' s'/ &
1378            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
1379            '                            bottom:     ', A, ' top:         ', A/&
1380            '       Maximum particle age:                 ',F9.1,' s'/ &
1381            '       Advection stopped at t = ',F9.1,' s'/)
1382481 FORMAT ('       Particles have random start positions'/)
1383482 FORMAT ('       Particles are advected only horizontally'/)
1384483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
1385484 FORMAT ('            Number of tails of the total domain: ',I10/ &
1386            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
1387            '            Maximum age of the end of the tail:  ',F8.2,' s')
1388485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
1389486 FORMAT ('       Particle statistics are written on file'/)
1390487 FORMAT ('       Number of particle groups: ',I2/)
1391488 FORMAT ('       SGS velocity components are used for particle advection'/ &
1392            '          minimum timestep for advection: ', F7.5/)
1393489 FORMAT ('       Number of particles simultaneously released at each ', &
1394                    'point: ', I5/)
1395490 FORMAT ('       Particle group ',I2,':'/ &
1396            '          Particle radius: ',E10.3, 'm')
1397491 FORMAT ('          Particle inertia is activated'/ &
1398            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
1399492 FORMAT ('          Particles are advected only passively (no inertia)'/)
1400493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
1401            '                                         y:',F8.1,' - ',F8.1,' m'/&
1402            '                                         z:',F8.1,' - ',F8.1,' m'/&
1403            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
1404                       ' m  dz = ',F8.1,' m'/)
1405494 FORMAT ('       Output of particle time series in NetCDF format every ', &
1406                    F8.2,' s'/)
1407495 FORMAT ('       Number of particles in total domain: ',I10/)
1408#endif
1409500 FORMAT (//' 1D-Model parameters:'/                           &
1410              ' -------------------'//                           &
1411            '    Simulation time:                   ',F8.1,' s'/ &
1412            '    Run-controll output every:         ',F8.1,' s'/ &
1413            '    Vertical profile output every:     ',F8.1,' s'/ &
1414            '    Mixing length calculation:         ',A/         &
1415            '    Dissipation calculation:           ',A/)
1416502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
1417
1418
1419 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.