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

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

preliminary update for changes concerning non-cyclic boundary conditions

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