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

Last change on this file since 419 was 411, checked in by heinze, 14 years ago

Large scale vertical motion (subsidence/ascent) can be applied to the prognostic equation for the potential temperature

  • Property svn:keywords set to Id
File size: 75.2 KB
RevLine 
[1]1 SUBROUTINE header
2
3!------------------------------------------------------------------------------!
[254]4! Current revisions:
[1]5! -----------------
[411]6! +large scale vertical motion (subsidence/ascent)
[392]7!
8!
[410]9! Branch revisions:
10! -----------------
11! Masked data output: + dt_domask, mask_01~20_x|y|z, mask_01~20_x|y|z_loop,
12! mask_scale|_x|y|z, masks, netcdf_format_mask[_av], skip_time_domask
13!
[392]14! Former revisions:
15! -----------------
16! $Id: header.f90 411 2009-12-11 14:15:58Z heinze $
17!
18! 346 2009-07-06 10:13:41Z raasch
[328]19! initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill'
[291]20! Coupling with independent precursor runs.
[254]21! Output of messages replaced by message handling routine.
[336]22! Output of several additional dvr parameters
[240]23! +canyon_height, canyon_width_x, canyon_width_y, canyon_wall_left,
[241]24! canyon_wall_south, conserve_volume_flow_mode, dp_external, dp_level_b,
25! dp_smooth, dpdxy, u_bulk, v_bulk
[256]26! topography_grid_convention moved from user_header
[292]27! small bugfix concerning 3d 64bit netcdf output format
[1]28!
[226]29! 206 2008-10-13 14:59:11Z raasch
30! Bugfix: error in zu index in case of section_xy = -1
31!
[200]32! 198 2008-09-17 08:55:28Z raasch
33! Format adjustments allowing output of larger revision numbers
34!
[198]35! 197 2008-09-16 15:29:03Z raasch
36! allow 100 spectra levels instead of 10 for consistency with
37! define_netcdf_header,
38! bugfix in the output of the characteristic levels of potential temperature,
39! geostrophic wind, scalar concentration, humidity and leaf area density,
40! output of turbulence recycling informations
41!
[139]42! 138 2007-11-28 10:03:58Z letzel
43! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
44! Allow two instead of one digit to specify isosurface and slicer variables.
45! Output of sorting frequency of particles
46!
[110]47! 108 2007-08-24 15:10:38Z letzel
48! Output of informations for coupled model runs (boundary conditions etc.)
49! + output of momentumfluxes at the top boundary
50! Rayleigh damping for ocean, e_init
51!
[98]52! 97 2007-06-21 08:23:15Z raasch
53! Adjustments for the ocean version.
54! use_pt_reference renamed use_reference
55!
[90]56! 87 2007-05-22 15:46:47Z raasch
57! Bugfix: output of use_upstream_for_tke
58!
[83]59! 82 2007-04-16 15:40:52Z raasch
60! Preprocessor strings for different linux clusters changed to "lc",
61! routine local_flush is used for buffer flushing
62!
[77]63! 76 2007-03-29 00:58:32Z raasch
64! Output of netcdf_64bit_3d, particles-package is now part of the default code,
65! output of the loop optimization method, moisture renamed humidity,
66! output of subversion revision number
67!
[39]68! 19 2007-02-23 04:53:48Z raasch
69! Output of scalar flux applied at top boundary
70!
[3]71! RCS Log replace by Id keyword, revision history cleaned up
72!
[1]73! Revision 1.63  2006/08/22 13:53:13  raasch
74! Output of dz_max
75!
76! Revision 1.1  1997/08/11 06:17:20  raasch
77! Initial revision
78!
79!
80! Description:
81! ------------
82! Writing a header with all important informations about the actual run.
83! This subroutine is called three times, two times at the beginning
84! (writing information on files RUN_CONTROL and HEADER) and one time at the
85! end of the run, then writing additional information about CPU-usage on file
86! header.
[411]87!-----------------------------------------------------------------------------!
[1]88
89    USE arrays_3d
90    USE control_parameters
91    USE cloud_parameters
92    USE cpulog
93    USE dvrp_variables
94    USE grid_variables
95    USE indices
96    USE model_1d
97    USE particle_attributes
98    USE pegrid
[411]99    USE subsidence_mod
[1]100    USE spectrum
101
102    IMPLICIT NONE
103
104    CHARACTER (LEN=1)  ::  prec
105    CHARACTER (LEN=2)  ::  do2d_mode
106    CHARACTER (LEN=5)  ::  section_chr
107    CHARACTER (LEN=9)  ::  time_to_string
108    CHARACTER (LEN=10) ::  coor_chr, host_chr
109    CHARACTER (LEN=16) ::  begin_chr
[200]110    CHARACTER (LEN=23) ::  ver_rev
[1]111    CHARACTER (LEN=40) ::  output_format
[167]112    CHARACTER (LEN=70) ::  char1, char2, dopr_chr, &
[1]113                           do2d_xy, do2d_xz, do2d_yz, do3d_chr, &
[410]114                           domask_chr, run_classification
[167]115    CHARACTER (LEN=86) ::  coordinates, gradients, learde, slices,  &
116                           temperatures, ugcomponent, vgcomponent
[1]117    CHARACTER (LEN=85) ::  roben, runten
118
[410]119    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)
120
121    INTEGER ::  av, bh, blx, bly, bxl, bxr, byn, bys, ch, count, cwx, cwy,  &
122         cxl, cxr, cyn, cys, dim, i, ihost, io, j, l, ll, m, mpi_type
[1]123    REAL    ::  cpuseconds_per_simulated_second
124
125!
126!-- Open the output file. At the end of the simulation, output is directed
127!-- to unit 19.
128    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
129         .NOT. simulated_time_at_begin /= simulated_time )  THEN
130       io = 15   !  header output on file RUN_CONTROL
131    ELSE
132       io = 19   !  header output on file HEADER
133    ENDIF
134    CALL check_open( io )
135
136!
137!-- At the end of the run, output file (HEADER) will be rewritten with
138!-- new informations
139    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
140
141!
142!-- Determine kind of model run
143    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
144       run_classification = '3D - restart run'
[328]145    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
146       run_classification = '3D - run with cyclic fill of 3D - prerun data'
[147]147    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
148       run_classification = '3D - run without 1D - prerun'
[197]149    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
[147]150       run_classification = '3D - run with 1D - prerun'
[197]151    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
152       run_classification = '3D - run initialized by user'
[1]153    ELSE
[254]154       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
155       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
[1]156    ENDIF
[97]157    IF ( ocean )  THEN
158       run_classification = 'ocean - ' // run_classification
159    ELSE
160       run_classification = 'atmosphere - ' // run_classification
161    ENDIF
[1]162
163!
164!-- Run-identification, date, time, host
165    host_chr = host(1:10)
[75]166    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
[102]167    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
[291]168    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
169#if defined( __mpi2 )
170       mpi_type = 2
171#else
172       mpi_type = 1
173#endif
174       WRITE ( io, 101 )  mpi_type, coupling_mode
175    ENDIF
[102]176    WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr, &
177                       ADJUSTR( host_chr )
[1]178#if defined( __parallel )
179    IF ( npex == -1  .AND.  pdims(2) /= 1 )  THEN
180       char1 = 'calculated'
181    ELSEIF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.  &
182               host(1:2) == 'lc' )  .AND.                          &
183             npex == -1  .AND.  pdims(2) == 1 )  THEN
184       char1 = 'forced'
185    ELSE
186       char1 = 'predefined'
187    ENDIF
188    IF ( threads_per_task == 1 )  THEN
[102]189       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
[1]190    ELSE
[102]191       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
[1]192                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
193    ENDIF
194    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
195           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
196         npex == -1  .AND.  pdims(2) == 1 )                      &
197    THEN
[102]198       WRITE ( io, 106 )
[1]199    ELSEIF ( pdims(2) == 1 )  THEN
[102]200       WRITE ( io, 107 )  'x'
[1]201    ELSEIF ( pdims(1) == 1 )  THEN
[102]202       WRITE ( io, 107 )  'y'
[1]203    ENDIF
[102]204    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
[1]205#endif
206    WRITE ( io, 99 )
207
208!
209!-- Numerical schemes
210    WRITE ( io, 110 )
211    IF ( psolver(1:7) == 'poisfft' )  THEN
212       WRITE ( io, 111 )  TRIM( fft_method )
213       IF ( psolver == 'poisfft_hybrid' )  WRITE ( io, 138 )
214    ELSEIF ( psolver == 'sor' )  THEN
215       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
216    ELSEIF ( psolver == 'multigrid' )  THEN
217       WRITE ( io, 135 )  cycle_mg, maximum_grid_level, ngsrb
218       IF ( mg_cycles == -1 )  THEN
219          WRITE ( io, 140 )  residual_limit
220       ELSE
221          WRITE ( io, 141 )  mg_cycles
222       ENDIF
223       IF ( mg_switch_to_pe0_level == 0 )  THEN
224          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
225                             nzt_mg(1)
[197]226       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
[1]227          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
228                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
229                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
230                             nzt_mg(mg_switch_to_pe0_level),    &
231                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
232                             nzt_mg(1)
233       ENDIF
234    ENDIF
235    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
236    THEN
237       WRITE ( io, 142 )
238    ENDIF
239
240    IF ( momentum_advec == 'pw-scheme' )  THEN
241       WRITE ( io, 113 )
242    ELSE
243       WRITE ( io, 114 )
244       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
245       IF ( overshoot_limit_u /= 0.0  .OR.  overshoot_limit_v /= 0.0  .OR. &
246            overshoot_limit_w /= 0.0 )  THEN
247          WRITE ( io, 127 )  overshoot_limit_u, overshoot_limit_v, &
248                             overshoot_limit_w
249       ENDIF
250       IF ( ups_limit_u /= 0.0  .OR.  ups_limit_v /= 0.0  .OR. &
251            ups_limit_w /= 0.0 )                               &
252       THEN
253          WRITE ( io, 125 )  ups_limit_u, ups_limit_v, ups_limit_w
254       ENDIF
255       IF ( long_filter_factor /= 0.0 )  WRITE ( io, 115 )  long_filter_factor
256    ENDIF
257    IF ( scalar_advec == 'pw-scheme' )  THEN
258       WRITE ( io, 116 )
259    ELSEIF ( scalar_advec == 'ups-scheme' )  THEN
260       WRITE ( io, 117 )
261       IF ( cut_spline_overshoot )  WRITE ( io, 124 )
262       IF ( overshoot_limit_e /= 0.0  .OR.  overshoot_limit_pt /= 0.0 )  THEN
263          WRITE ( io, 128 )  overshoot_limit_e, overshoot_limit_pt
264       ENDIF
265       IF ( ups_limit_e /= 0.0  .OR.  ups_limit_pt /= 0.0 )  THEN
266          WRITE ( io, 126 )  ups_limit_e, ups_limit_pt
267       ENDIF
268    ELSE
269       WRITE ( io, 118 )
270    ENDIF
[63]271
272    WRITE ( io, 139 )  TRIM( loop_optimization )
273
[1]274    IF ( galilei_transformation )  THEN
275       IF ( use_ug_for_galilei_tr )  THEN
276          char1 = 'geostrophic wind'
277       ELSE
278          char1 = 'mean wind in model domain'
279       ENDIF
280       IF ( simulated_time_at_begin == simulated_time )  THEN
281          char2 = 'at the start of the run'
282       ELSE
283          char2 = 'at the end of the run'
284       ENDIF
285       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ), &
286                          advected_distance_x/1000.0, advected_distance_y/1000.0
287    ENDIF
288    IF ( timestep_scheme == 'leapfrog' )  THEN
289       WRITE ( io, 120 )
290    ELSEIF ( timestep_scheme == 'leapfrog+euler' )  THEN
291       WRITE ( io, 121 )
292    ELSE
293       WRITE ( io, 122 )  timestep_scheme
294    ENDIF
[87]295    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
[1]296    IF ( rayleigh_damping_factor /= 0.0 )  THEN
[108]297       IF ( .NOT. ocean )  THEN
298          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
299               rayleigh_damping_factor
300       ELSE
301          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
302               rayleigh_damping_factor
303       ENDIF
[1]304    ENDIF
[75]305    IF ( humidity )  THEN
[1]306       IF ( .NOT. cloud_physics )  THEN
307          WRITE ( io, 129 )
308       ELSE
309          WRITE ( io, 130 )
310          WRITE ( io, 131 )
311          IF ( radiation )      WRITE ( io, 132 )
312          IF ( precipitation )  WRITE ( io, 133 )
313       ENDIF
314    ENDIF
315    IF ( passive_scalar )  WRITE ( io, 134 )
[240]316    IF ( conserve_volume_flow )  THEN
[241]317       WRITE ( io, 150 )  conserve_volume_flow_mode
318       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
319          WRITE ( io, 151 )  u_bulk, v_bulk
320       ENDIF
[240]321    ELSEIF ( dp_external )  THEN
322       IF ( dp_smooth )  THEN
[241]323          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
[240]324       ELSE
[241]325          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
[240]326       ENDIF
327    ENDIF
[411]328    IF ( large_scale_subsidence )  THEN
329        WRITE ( io, 153 )
330        WRITE ( io, 154 )
331    ENDIF
[1]332    WRITE ( io, 99 )
333
334!
335!-- Runtime and timestep informations
336    WRITE ( io, 200 )
337    IF ( .NOT. dt_fixed )  THEN
338       WRITE ( io, 201 )  dt_max, cfl_factor
339    ELSE
340       WRITE ( io, 202 )  dt
341    ENDIF
342    WRITE ( io, 203 )  simulated_time_at_begin, end_time
343
344    IF ( time_restart /= 9999999.9  .AND. &
345         simulated_time_at_begin == simulated_time )  THEN
346       IF ( dt_restart == 9999999.9 )  THEN
347          WRITE ( io, 204 )  ' Restart at:       ',time_restart
348       ELSE
349          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
350       ENDIF
351    ENDIF
352
353    IF ( simulated_time_at_begin /= simulated_time )  THEN
354       i = MAX ( log_point_s(10)%counts, 1 )
355       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0 )  THEN
356          cpuseconds_per_simulated_second = 0.0
357       ELSE
358          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
359                                            ( simulated_time -    &
360                                              simulated_time_at_begin )
361       ENDIF
362       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum, &
363                          log_point_s(10)%sum / REAL( i ),     &
364                          cpuseconds_per_simulated_second
365       IF ( time_restart /= 9999999.9  .AND.  time_restart < end_time )  THEN
366          IF ( dt_restart == 9999999.9 )  THEN
367             WRITE ( io, 204 )  ' Next restart at:  ',time_restart
368          ELSE
369             WRITE ( io, 205 )  ' Next restart at:  ',time_restart, dt_restart
370          ENDIF
371       ENDIF
372    ENDIF
373
374!
[291]375!-- Start time for coupled runs, if independent precursor runs for atmosphere
376!-- and ocean are used. In this case, coupling_start_time defines the time
377!-- when the coupling is switched on.
378    IF ( coupling_start_time /= 0.0 )  THEN
379       IF ( coupling_start_time >= simulated_time_at_begin )  THEN
380          char1 = 'Precursor run for a coupled atmosphere-ocean run'
381       ELSE
382          char1 = 'Coupled atmosphere-ocean run following independent ' // &
383                  'precursor runs'
384       ENDIF
385       WRITE ( io, 207 )  char1, coupling_start_time
386    ENDIF
387
388!
[1]389!-- Computational grid
[94]390    IF ( .NOT. ocean )  THEN
391       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
392       IF ( dz_stretch_level_index < nzt+1 )  THEN
393          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
394                             dz_stretch_factor, dz_max
395       ENDIF
396    ELSE
397       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
398       IF ( dz_stretch_level_index > 0 )  THEN
399          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
400                             dz_stretch_factor, dz_max
401       ENDIF
[1]402    ENDIF
403    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
404                       MIN( nnz+2, nzt+2 )
[76]405    IF ( numprocs > 1 )  THEN
406       IF ( nxa == nx  .AND.  nya == ny  .AND.  nza == nz )  THEN
407          WRITE ( io, 255 )
408       ELSE
409          WRITE ( io, 256 )  nnx-(nxa-nx), nny-(nya-ny), nzt+2
410       ENDIF
[1]411    ENDIF
412    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
413
414!
415!-- Topography
416    WRITE ( io, 270 )  topography
417    SELECT CASE ( TRIM( topography ) )
418
419       CASE ( 'flat' )
420          ! no actions necessary
421
422       CASE ( 'single_building' )
423          blx = INT( building_length_x / dx )
424          bly = INT( building_length_y / dy )
425          bh  = INT( building_height / dz )
426
427          IF ( building_wall_left == 9999999.9 )  THEN
428             building_wall_left = ( nx + 1 - blx ) / 2 * dx
429          ENDIF
430          bxl = INT ( building_wall_left / dx + 0.5 )
431          bxr = bxl + blx
432
433          IF ( building_wall_south == 9999999.9 )  THEN
434             building_wall_south = ( ny + 1 - bly ) / 2 * dy
435          ENDIF
436          bys = INT ( building_wall_south / dy + 0.5 )
437          byn = bys + bly
438
439          WRITE ( io, 271 )  building_length_x, building_length_y, &
440                             building_height, bxl, bxr, bys, byn
441
[240]442       CASE ( 'single_street_canyon' )
443          ch  = NINT( canyon_height / dz )
444          IF ( canyon_width_x /= 9999999.9 )  THEN
445!
446!--          Street canyon in y direction
447             cwx = NINT( canyon_width_x / dx )
448             IF ( canyon_wall_left == 9999999.9 )  THEN
449                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
450             ENDIF
451             cxl = NINT( canyon_wall_left / dx )
452             cxr = cxl + cwx
453             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
454
455          ELSEIF ( canyon_width_y /= 9999999.9 )  THEN
456!
457!--          Street canyon in x direction
458             cwy = NINT( canyon_width_y / dy )
459             IF ( canyon_wall_south == 9999999.9 )  THEN
460                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
461             ENDIF
462             cys = NINT( canyon_wall_south / dy )
463             cyn = cys + cwy
464             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
465          ENDIF
466
[1]467    END SELECT
468
[256]469    IF ( TRIM( topography ) /= 'flat' )  THEN
470       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
471          IF ( TRIM( topography ) == 'single_building' .OR.  &
472               TRIM( topography ) == 'single_street_canyon' )  THEN
473             WRITE ( io, 278 )
474          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
475             WRITE ( io, 279 )
476          ENDIF
477       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
478          WRITE ( io, 278 )
479       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
480          WRITE ( io, 279 )
481       ENDIF
482    ENDIF
483
[138]484    IF ( plant_canopy ) THEN
485
486       WRITE ( io, 280 ) canopy_mode, pch_index, drag_coefficient
[153]487       IF ( passive_scalar ) THEN
488          WRITE ( io, 281 ) scalar_exchange_coefficient,   &
489                            leaf_surface_concentration
490       ENDIF
[138]491
[1]492!
[153]493!--    Heat flux at the top of vegetation
494       WRITE ( io, 282 ) cthf
495
496!
[138]497!--    Leaf area density profile
498!--    Building output strings, starting with surface value
499       WRITE ( learde, '(F6.2)' )  lad_surface
500       gradients = '------'
501       slices = '     0'
502       coordinates = '   0.0'
503       i = 1
504       DO  WHILE ( lad_vertical_gradient_level_ind(i) /= -9999 )
505
506          WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
507          learde = TRIM( learde ) // ' ' // TRIM( coor_chr )
508
509          WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
510          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
511
512          WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
513          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
514
515          WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
516          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
517
518          i = i + 1
519       ENDDO
520
[153]521       WRITE ( io, 283 )  TRIM( coordinates ), TRIM( learde ), &
[138]522                          TRIM( gradients ), TRIM( slices )
523
524    ENDIF
525
526!
[1]527!-- Boundary conditions
528    IF ( ibc_p_b == 0 )  THEN
529       runten = 'p(0)     = 0      |'
530    ELSEIF ( ibc_p_b == 1 )  THEN
531       runten = 'p(0)     = p(1)   |'
532    ELSE
533       runten = 'p(0)     = p(1) +R|'
534    ENDIF
535    IF ( ibc_p_t == 0 )  THEN
536       roben  = 'p(nzt+1) = 0      |'
537    ELSE
538       roben  = 'p(nzt+1) = p(nzt) |'
539    ENDIF
540
541    IF ( ibc_uv_b == 0 )  THEN
542       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
543    ELSE
544       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
545    ENDIF
[132]546    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
547       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
548    ELSEIF ( ibc_uv_t == 0 )  THEN
[1]549       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
550    ELSE
551       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
552    ENDIF
553
554    IF ( ibc_pt_b == 0 )  THEN
555       runten = TRIM( runten ) // ' pt(0)   = pt_surface'
[102]556    ELSEIF ( ibc_pt_b == 1 )  THEN
[1]557       runten = TRIM( runten ) // ' pt(0)   = pt(1)'
[102]558    ELSEIF ( ibc_pt_b == 2 )  THEN
559       runten = TRIM( runten ) // ' pt(0) = from coupled model'
[1]560    ENDIF
561    IF ( ibc_pt_t == 0 )  THEN
[19]562       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
563    ELSEIF( ibc_pt_t == 1 )  THEN
564       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
565    ELSEIF( ibc_pt_t == 2 )  THEN
566       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
[1]567    ENDIF
568
569    WRITE ( io, 300 )  runten, roben
570
571    IF ( .NOT. constant_diffusion )  THEN
572       IF ( ibc_e_b == 1 )  THEN
573          runten = 'e(0)     = e(1)'
574       ELSE
575          runten = 'e(0)     = e(1) = (u*/0.1)**2'
576       ENDIF
577       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
578
[97]579       WRITE ( io, 301 )  'e', runten, roben       
[1]580
581    ENDIF
582
[97]583    IF ( ocean )  THEN
584       runten = 'sa(0)    = sa(1)'
585       IF ( ibc_sa_t == 0 )  THEN
586          roben =  'sa(nzt+1) = sa_surface'
[1]587       ELSE
[97]588          roben =  'sa(nzt+1) = sa(nzt)'
[1]589       ENDIF
[97]590       WRITE ( io, 301 ) 'sa', runten, roben
591    ENDIF
[1]592
[97]593    IF ( humidity )  THEN
594       IF ( ibc_q_b == 0 )  THEN
595          runten = 'q(0)     = q_surface'
596       ELSE
597          runten = 'q(0)     = q(1)'
598       ENDIF
599       IF ( ibc_q_t == 0 )  THEN
600          roben =  'q(nzt)   = q_top'
601       ELSE
602          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
603       ENDIF
604       WRITE ( io, 301 ) 'q', runten, roben
605    ENDIF
[1]606
[97]607    IF ( passive_scalar )  THEN
608       IF ( ibc_q_b == 0 )  THEN
609          runten = 's(0)     = s_surface'
610       ELSE
611          runten = 's(0)     = s(1)'
612       ENDIF
613       IF ( ibc_q_t == 0 )  THEN
614          roben =  's(nzt)   = s_top'
615       ELSE
616          roben =  's(nzt)   = s(nzt-1) + ds/dz'
617       ENDIF
618       WRITE ( io, 301 ) 's', runten, roben
[1]619    ENDIF
620
621    IF ( use_surface_fluxes )  THEN
622       WRITE ( io, 303 )
623       IF ( constant_heatflux )  THEN
624          WRITE ( io, 306 )  surface_heatflux
625          IF ( random_heatflux )  WRITE ( io, 307 )
626       ENDIF
[75]627       IF ( humidity  .AND.  constant_waterflux )  THEN
[1]628          WRITE ( io, 311 ) surface_waterflux
629       ENDIF
630       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
631          WRITE ( io, 313 ) surface_waterflux
632       ENDIF
633    ENDIF
634
[19]635    IF ( use_top_fluxes )  THEN
636       WRITE ( io, 304 )
[102]637       IF ( coupling_mode == 'uncoupled' )  THEN
[151]638          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
[102]639          IF ( constant_top_heatflux )  THEN
640             WRITE ( io, 306 )  top_heatflux
641          ENDIF
642       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
643          WRITE ( io, 316 )
[19]644       ENDIF
[97]645       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
646          WRITE ( io, 309 )  top_salinityflux
647       ENDIF
[75]648       IF ( humidity  .OR.  passive_scalar )  THEN
[19]649          WRITE ( io, 315 )
650       ENDIF
651    ENDIF
652
[1]653    IF ( prandtl_layer )  THEN
[94]654       WRITE ( io, 305 )  0.5 * (zu(1)-zu(0)), roughness_length, kappa, &
655                          rif_min, rif_max
[1]656       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
[75]657       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
[1]658          WRITE ( io, 312 )
659       ENDIF
660       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
661          WRITE ( io, 314 )
662       ENDIF
663    ELSE
664       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
665          WRITE ( io, 310 )  rif_min, rif_max
666       ENDIF
667    ENDIF
668
669    WRITE ( io, 317 )  bc_lr, bc_ns
670    IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
671       WRITE ( io, 318 )  outflow_damping_width, km_damp_max
[151]672       IF ( turbulent_inflow )  THEN
673          WRITE ( io, 319 )  recycling_width, recycling_plane, &
674                             inflow_damping_height, inflow_damping_width
675       ENDIF
[1]676    ENDIF
677
678!
679!-- Listing of 1D-profiles
[151]680    WRITE ( io, 325 )  dt_dopr_listing
[1]681    IF ( averaging_interval_pr /= 0.0 )  THEN
[151]682       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]683    ENDIF
684
685!
686!-- DATA output
687    WRITE ( io, 330 )
688    IF ( averaging_interval_pr /= 0.0 )  THEN
[151]689       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]690    ENDIF
691
692!
693!-- 1D-profiles
[346]694    dopr_chr = 'Profile:'
[1]695    IF ( dopr_n /= 0 )  THEN
696       WRITE ( io, 331 )
697
698       output_format = ''
699       IF ( netcdf_output )  THEN
700          IF ( netcdf_64bit )  THEN
701             output_format = 'netcdf (64 bit offset)'
702          ELSE
703             output_format = 'netcdf'
704          ENDIF
705       ENDIF
706       IF ( profil_output )  THEN
707          IF ( netcdf_output )  THEN
708             output_format = TRIM( output_format ) // ' and profil'
709          ELSE
710             output_format = 'profil'
711          ENDIF
712       ENDIF
[292]713       WRITE ( io, 344 )  output_format
[1]714
715       DO  i = 1, dopr_n
716          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
717          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
718             WRITE ( io, 332 )  dopr_chr
719             dopr_chr = '       :'
720          ENDIF
721       ENDDO
722
723       IF ( dopr_chr /= '' )  THEN
724          WRITE ( io, 332 )  dopr_chr
725       ENDIF
726       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
727       IF ( skip_time_dopr /= 0.0 )  WRITE ( io, 339 )  skip_time_dopr
728    ENDIF
729
730!
731!-- 2D-arrays
732    DO  av = 0, 1
733
734       i = 1
735       do2d_xy = ''
736       do2d_xz = ''
737       do2d_yz = ''
738       DO  WHILE ( do2d(av,i) /= ' ' )
739
740          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
741          do2d_mode = do2d(av,i)(l-1:l)
742
743          SELECT CASE ( do2d_mode )
744             CASE ( 'xy' )
745                ll = LEN_TRIM( do2d_xy )
746                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
747             CASE ( 'xz' )
748                ll = LEN_TRIM( do2d_xz )
749                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
750             CASE ( 'yz' )
751                ll = LEN_TRIM( do2d_yz )
752                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
753          END SELECT
754
755          i = i + 1
756
757       ENDDO
758
759       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
760              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
761              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) )  .AND. &
762            ( netcdf_output  .OR.  iso2d_output ) )  THEN
763
764          IF (  av == 0 )  THEN
765             WRITE ( io, 334 )  ''
766          ELSE
767             WRITE ( io, 334 )  '(time-averaged)'
768          ENDIF
769
770          IF ( do2d_at_begin )  THEN
771             begin_chr = 'and at the start'
772          ELSE
773             begin_chr = ''
774          ENDIF
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 ( iso2d_output )  THEN
785             IF ( netcdf_output )  THEN
786                output_format = TRIM( output_format ) // ' and iso2d'
787             ELSE
788                output_format = 'iso2d'
789             ENDIF
790          ENDIF
[292]791          WRITE ( io, 344 )  output_format
[1]792
793          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
794             i = 1
795             slices = '/'
796             coordinates = '/'
797!
798!--          Building strings with index and coordinate informations of the
799!--          slices
800             DO  WHILE ( section(i,1) /= -9999 )
801
802                WRITE (section_chr,'(I5)')  section(i,1)
803                section_chr = ADJUSTL( section_chr )
804                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
805
[206]806                IF ( section(i,1) == -1 )  THEN
807                   WRITE (coor_chr,'(F10.1)')  -1.0
808                ELSE
809                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
810                ENDIF
[1]811                coor_chr = ADJUSTL( coor_chr )
812                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
813
814                i = i + 1
815             ENDDO
816             IF ( av == 0 )  THEN
817                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
818                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
819                                   TRIM( coordinates )
820                IF ( skip_time_do2d_xy /= 0.0 )  THEN
821                   WRITE ( io, 339 )  skip_time_do2d_xy
822                ENDIF
823             ELSE
824                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
825                                   TRIM( begin_chr ), averaging_interval, &
826                                   dt_averaging_input, 'k', TRIM( slices ), &
827                                   TRIM( coordinates )
828                IF ( skip_time_data_output_av /= 0.0 )  THEN
829                   WRITE ( io, 339 )  skip_time_data_output_av
830                ENDIF
831             ENDIF
832
833          ENDIF
834
835          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
836             i = 1
837             slices = '/'
838             coordinates = '/'
839!
840!--          Building strings with index and coordinate informations of the
841!--          slices
842             DO  WHILE ( section(i,2) /= -9999 )
843
844                WRITE (section_chr,'(I5)')  section(i,2)
845                section_chr = ADJUSTL( section_chr )
846                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
847
848                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
849                coor_chr = ADJUSTL( coor_chr )
850                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
851
852                i = i + 1
853             ENDDO
854             IF ( av == 0 )  THEN
855                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
856                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
857                                   TRIM( coordinates )
858                IF ( skip_time_do2d_xz /= 0.0 )  THEN
859                   WRITE ( io, 339 )  skip_time_do2d_xz
860                ENDIF
861             ELSE
862                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
863                                   TRIM( begin_chr ), averaging_interval, &
864                                   dt_averaging_input, 'j', TRIM( slices ), &
865                                   TRIM( coordinates )
866                IF ( skip_time_data_output_av /= 0.0 )  THEN
867                   WRITE ( io, 339 )  skip_time_data_output_av
868                ENDIF
869             ENDIF
870          ENDIF
871
872          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
873             i = 1
874             slices = '/'
875             coordinates = '/'
876!
877!--          Building strings with index and coordinate informations of the
878!--          slices
879             DO  WHILE ( section(i,3) /= -9999 )
880
881                WRITE (section_chr,'(I5)')  section(i,3)
882                section_chr = ADJUSTL( section_chr )
883                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
884
885                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
886                coor_chr = ADJUSTL( coor_chr )
887                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
888
889                i = i + 1
890             ENDDO
891             IF ( av == 0 )  THEN
892                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
893                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
894                                   TRIM( coordinates )
895                IF ( skip_time_do2d_yz /= 0.0 )  THEN
896                   WRITE ( io, 339 )  skip_time_do2d_yz
897                ENDIF
898             ELSE
899                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
900                                   TRIM( begin_chr ), averaging_interval, &
901                                   dt_averaging_input, 'i', TRIM( slices ), &
902                                   TRIM( coordinates )
903                IF ( skip_time_data_output_av /= 0.0 )  THEN
904                   WRITE ( io, 339 )  skip_time_data_output_av
905                ENDIF
906             ENDIF
907          ENDIF
908
909       ENDIF
910
911    ENDDO
912
913!
914!-- 3d-arrays
915    DO  av = 0, 1
916
917       i = 1
918       do3d_chr = ''
919       DO  WHILE ( do3d(av,i) /= ' ' )
920
921          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
922          i = i + 1
923
924       ENDDO
925
926       IF ( do3d_chr /= '' )  THEN
927          IF ( av == 0 )  THEN
928             WRITE ( io, 336 )  ''
929          ELSE
930             WRITE ( io, 336 )  '(time-averaged)'
931          ENDIF
932
933          output_format = ''
934          IF ( netcdf_output )  THEN
[292]935             IF ( netcdf_64bit_3d )  THEN
[1]936                output_format = 'netcdf (64 bit offset)'
937             ELSE
938                output_format = 'netcdf'
939             ENDIF
940          ENDIF
941          IF ( avs_output )  THEN
942             IF ( netcdf_output )  THEN
943                output_format = TRIM( output_format ) // ' and avs'
944             ELSE
945                output_format = 'avs'
946             ENDIF
947          ENDIF
[292]948          WRITE ( io, 344 )  output_format
[1]949
950          IF ( do3d_at_begin )  THEN
951             begin_chr = 'and at the start'
952          ELSE
953             begin_chr = ''
954          ENDIF
955          IF ( av == 0 )  THEN
956             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
957                                zu(nz_do3d), nz_do3d
958          ELSE
959             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
960                                TRIM( begin_chr ), averaging_interval, &
961                                dt_averaging_input, zu(nz_do3d), nz_do3d
962          ENDIF
963
964          IF ( do3d_compress )  THEN
965             do3d_chr = ''
966             i = 1
967             DO WHILE ( do3d(av,i) /= ' ' )
968
969                SELECT CASE ( do3d(av,i) )
970                   CASE ( 'u' )
971                      j = 1
972                   CASE ( 'v' )
973                      j = 2
974                   CASE ( 'w' )
975                      j = 3
976                   CASE ( 'p' )
977                      j = 4
978                   CASE ( 'pt' )
979                      j = 5
980                END SELECT
981                WRITE ( prec, '(I1)' )  plot_3d_precision(j)%precision
982                do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // &
983                           ':' // prec // ','
984                i = i + 1
985
986             ENDDO
987             WRITE ( io, 338 )  do3d_chr
988
989          ENDIF
990
991          IF ( av == 0 )  THEN
992             IF ( skip_time_do3d /= 0.0 )  THEN
993                WRITE ( io, 339 )  skip_time_do3d
994             ENDIF
995          ELSE
996             IF ( skip_time_data_output_av /= 0.0 )  THEN
997                WRITE ( io, 339 )  skip_time_data_output_av
998             ENDIF
999          ENDIF
1000
1001       ENDIF
1002
1003    ENDDO
1004
1005!
[410]1006!-- masked arrays
1007    IF ( masks > 0 )  WRITE ( io, 345 )  &
1008         mask_scale_x, mask_scale_y, mask_scale_z
1009    DO  mid = 1, masks
1010       DO  av = 0, 1
1011
1012          i = 1
1013          domask_chr = ''
1014          DO  WHILE ( domask(mid,av,i) /= ' ' )
1015             domask_chr = TRIM( domask_chr ) // ' ' //  &
1016                          TRIM( domask(mid,av,i) ) // ','
1017             i = i + 1
1018          ENDDO
1019
1020          IF ( domask_chr /= '' )  THEN
1021             IF ( av == 0 )  THEN
1022                WRITE ( io, 346 )  '', mid
1023             ELSE
1024                WRITE ( io, 346 )  ' (time-averaged)', mid
1025             ENDIF
1026
1027             output_format = ''
1028             IF ( netcdf_output )  THEN
1029                SELECT CASE ( nc_format_mask(mid,av) )
1030                   CASE ( 1 )
1031                      output_format = 'netcdf (classic format)'
1032                   CASE ( 2 )
1033                      output_format = 'netcdf (64bit offset format)'
1034                   CASE ( 3 )
1035                      output_format = 'netcdf (NetCDF 4 format)'
1036                   CASE ( 4 )
1037                      output_format = 'netcdf (NetCDF 4 classic model format)'
1038                END SELECT
1039             ENDIF
1040             WRITE ( io, 344 )  output_format
1041
1042             IF ( av == 0 )  THEN
1043                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1044             ELSE
1045                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1046                                   averaging_interval, dt_averaging_input
1047             ENDIF
1048
1049             IF ( av == 0 )  THEN
1050                IF ( skip_time_domask(mid) /= 0.0 )  THEN
1051                   WRITE ( io, 339 )  skip_time_domask(mid)
1052                ENDIF
1053             ELSE
1054                IF ( skip_time_data_output_av /= 0.0 )  THEN
1055                   WRITE ( io, 339 )  skip_time_data_output_av
1056                ENDIF
1057             ENDIF
1058!
1059!--          output locations
1060             DO  dim = 1, 3
1061                IF ( mask(mid,dim,1) >= 0.0 )  THEN
1062                   count = 0
1063                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0 )
1064                      count = count + 1
1065                   ENDDO
1066                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1067                                      mask(mid,dim,:count)
1068                ELSEIF ( mask_loop(mid,dim,1) < 0.0 .AND.  &
1069                         mask_loop(mid,dim,2) < 0.0 .AND.  &
1070                         mask_loop(mid,dim,3) == 0.0 )  THEN
1071                   WRITE ( io, 350 )  dir(dim), dir(dim)
1072                ELSEIF ( mask_loop(mid,dim,3) == 0.0 )  THEN
1073                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1074                                      mask_loop(mid,dim,1:2)
1075                ELSE
1076                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1077                                      mask_loop(mid,dim,1:3)
1078                ENDIF
1079             ENDDO
1080          ENDIF
1081
1082       ENDDO
1083    ENDDO
1084
1085!
[1]1086!-- Timeseries
1087    IF ( dt_dots /= 9999999.9 )  THEN
1088       WRITE ( io, 340 )
1089
1090       output_format = ''
1091       IF ( netcdf_output )  THEN
1092          IF ( netcdf_64bit )  THEN
1093             output_format = 'netcdf (64 bit offset)'
1094          ELSE
1095             output_format = 'netcdf'
1096          ENDIF
1097       ENDIF
1098       IF ( profil_output )  THEN
1099          IF ( netcdf_output )  THEN
1100             output_format = TRIM( output_format ) // ' and profil'
1101          ELSE
1102             output_format = 'profil'
1103          ENDIF
1104       ENDIF
[292]1105       WRITE ( io, 344 )  output_format
[1]1106       WRITE ( io, 341 )  dt_dots
1107    ENDIF
1108
1109#if defined( __dvrp_graphics )
1110!
1111!-- Dvrp-output
1112    IF ( dt_dvrp /= 9999999.9 )  THEN
1113       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1114                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1115       i = 1
1116       l = 0
[336]1117       m = 0
[1]1118       DO WHILE ( mode_dvrp(i) /= ' ' )
1119          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
[130]1120             READ ( mode_dvrp(i), '(10X,I2)' )  j
[1]1121             l = l + 1
1122             IF ( do3d(0,j) /= ' ' )  THEN
[336]1123                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1124                                   isosurface_color(:,l)
[1]1125             ENDIF
1126          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
[130]1127             READ ( mode_dvrp(i), '(6X,I2)' )  j
[336]1128             m = m + 1
1129             IF ( do2d(0,j) /= ' ' )  THEN
1130                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1131                                   slicer_range_limits_dvrp(:,m)
1132             ENDIF
[1]1133          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
[336]1134             WRITE ( io, 363 )  dvrp_psize
1135             IF ( particle_dvrpsize /= 'none' )  THEN
1136                WRITE ( io, 364 )  'size', TRIM( particle_dvrpsize ), &
1137                                   dvrpsize_interval
1138             ENDIF
1139             IF ( particle_color /= 'none' )  THEN
1140                WRITE ( io, 364 )  'color', TRIM( particle_color ), &
1141                                   color_interval
1142             ENDIF
[1]1143          ENDIF
1144          i = i + 1
1145       ENDDO
[237]1146
[336]1147       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1148                          superelevation_y, superelevation, clip_dvrp_l, &
1149                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1150
1151       IF ( TRIM( topography ) /= 'flat' )  THEN
1152          WRITE ( io, 366 )  topography_color
1153          IF ( cluster_size > 1 )  THEN
1154             WRITE ( io, 367 )  cluster_size
1155          ENDIF
[237]1156       ENDIF
1157
[1]1158    ENDIF
1159#endif
1160
1161#if defined( __spectra )
1162!
1163!-- Spectra output
1164    IF ( dt_dosp /= 9999999.9 ) THEN
1165       WRITE ( io, 370 )
1166
1167       output_format = ''
1168       IF ( netcdf_output )  THEN
1169          IF ( netcdf_64bit )  THEN
1170             output_format = 'netcdf (64 bit offset)'
1171          ELSE
1172             output_format = 'netcdf'
1173          ENDIF
1174       ENDIF
1175       IF ( profil_output )  THEN
1176          IF ( netcdf_output )  THEN
1177             output_format = TRIM( output_format ) // ' and profil'
1178          ELSE
1179             output_format = 'profil'
1180          ENDIF
1181       ENDIF
[292]1182       WRITE ( io, 344 )  output_format
[1]1183       WRITE ( io, 371 )  dt_dosp
1184       IF ( skip_time_dosp /= 0.0 )  WRITE ( io, 339 )  skip_time_dosp
1185       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
1186                          ( spectra_direction(i), i = 1,10 ),  &
[189]1187                          ( comp_spectra_level(i), i = 1,100 ), &
1188                          ( plot_spectra_level(i), i = 1,100 ), &
[1]1189                          averaging_interval_sp, dt_averaging_input_pr
1190    ENDIF
1191#endif
1192
1193    WRITE ( io, 99 )
1194
1195!
1196!-- Physical quantities
1197    WRITE ( io, 400 )
1198
1199!
1200!-- Geostrophic parameters
1201    WRITE ( io, 410 )  omega, phi, f, fs
1202
1203!
1204!-- Other quantities
1205    WRITE ( io, 411 )  g
[97]1206    IF ( use_reference )  THEN
1207       IF ( ocean )  THEN
1208          WRITE ( io, 412 )  prho_reference
1209       ELSE
1210          WRITE ( io, 413 )  pt_reference
1211       ENDIF
1212    ENDIF
[1]1213
1214!
1215!-- Cloud physics parameters
1216    IF ( cloud_physics ) THEN
[57]1217       WRITE ( io, 415 )
1218       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
[1]1219    ENDIF
1220
1221!-- Profile of the geostrophic wind (component ug)
1222!-- Building output strings
1223    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
1224    gradients = '------'
1225    slices = '     0'
1226    coordinates = '   0.0'
1227    i = 1
1228    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
1229     
[167]1230       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
[1]1231       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
1232
[167]1233       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
[1]1234       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1235
[167]1236       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
[1]1237       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1238
[167]1239       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
[1]1240       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1241
1242       i = i + 1
1243    ENDDO
1244
1245    WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
1246                       TRIM( gradients ), TRIM( slices )
1247
1248!-- Profile of the geostrophic wind (component vg)
1249!-- Building output strings
1250    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
1251    gradients = '------'
1252    slices = '     0'
1253    coordinates = '   0.0'
1254    i = 1
1255    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
1256
[167]1257       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
[1]1258       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
1259
[167]1260       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
[1]1261       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1262
[167]1263       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
[1]1264       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1265
[167]1266       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
[1]1267       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
1268
1269       i = i + 1 
1270    ENDDO
1271
1272    WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
1273                       TRIM( gradients ), TRIM( slices )
1274
1275!
1276!-- Initial temperature profile
1277!-- Building output strings, starting with surface temperature
1278    WRITE ( temperatures, '(F6.2)' )  pt_surface
1279    gradients = '------'
1280    slices = '     0'
1281    coordinates = '   0.0'
1282    i = 1
1283    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1284
[94]1285       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1286       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
[1]1287
[94]1288       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1289       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
[1]1290
[94]1291       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1292       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
[1]1293
[94]1294       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1295       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
[1]1296
1297       i = i + 1
1298    ENDDO
1299
1300    WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1301                       TRIM( gradients ), TRIM( slices )
1302
1303!
1304!-- Initial humidity profile
1305!-- Building output strings, starting with surface humidity
[75]1306    IF ( humidity  .OR.  passive_scalar )  THEN
[1]1307       WRITE ( temperatures, '(E8.1)' )  q_surface
1308       gradients = '--------'
1309       slices = '       0'
1310       coordinates = '     0.0'
1311       i = 1
1312       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1313         
1314          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1315          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1316
1317          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1318          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1319         
1320          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1321          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1322         
1323          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1324          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1325
1326          i = i + 1
1327       ENDDO
1328
[75]1329       IF ( humidity )  THEN
[1]1330          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1331                             TRIM( gradients ), TRIM( slices )
1332       ELSE
1333          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1334                             TRIM( gradients ), TRIM( slices )
1335       ENDIF
1336    ENDIF
1337
1338!
[97]1339!-- Initial salinity profile
1340!-- Building output strings, starting with surface salinity
1341    IF ( ocean )  THEN
1342       WRITE ( temperatures, '(F6.2)' )  sa_surface
1343       gradients = '------'
1344       slices = '     0'
1345       coordinates = '   0.0'
1346       i = 1
1347       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1348
1349          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1350          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1351
1352          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1353          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1354
1355          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1356          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1357
1358          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1359          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1360
1361          i = i + 1
1362       ENDDO
1363
1364       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1365                          TRIM( gradients ), TRIM( slices )
1366    ENDIF
1367
1368!
[411]1369!-- Profile for the large scale vertial velocity
1370!-- Building output strings, starting with surface value
1371    IF ( large_scale_subsidence )  THEN
1372       temperatures = '   0.0'
1373       gradients = '------'
1374       slices = '     0'
1375       coordinates = '   0.0'
1376       i = 1
1377       DO  WHILE ( ws_vertical_gradient_level_ind(i) /= -9999 )
1378
1379          WRITE (coor_chr,'(E10.2,7X)')  &
1380                                w_subs(ws_vertical_gradient_level_ind(i))
1381          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1382
1383          WRITE (coor_chr,'(E10.2,7X)')  ws_vertical_gradient(i)
1384          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1385
1386          WRITE (coor_chr,'(I10,7X)')  ws_vertical_gradient_level_ind(i)
1387          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1388
1389          WRITE (coor_chr,'(F10.2,7X)')  ws_vertical_gradient_level(i)
1390          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1391
1392          i = i + 1
1393       ENDDO
1394
1395       WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
1396                          TRIM( gradients ), TRIM( slices )
1397    ENDIF
1398
1399!
[1]1400!-- LES / turbulence parameters
1401    WRITE ( io, 450 )
1402
1403!--
1404! ... LES-constants used must still be added here
1405!--
1406    IF ( constant_diffusion )  THEN
1407       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1408                          prandtl_number
1409    ENDIF
1410    IF ( .NOT. constant_diffusion)  THEN
[108]1411       IF ( e_init > 0.0 )  WRITE ( io, 455 )  e_init
[1]1412       IF ( e_min > 0.0 )  WRITE ( io, 454 )  e_min
1413       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1414       IF ( adjust_mixing_length  .AND.  prandtl_layer )  WRITE ( io, 452 )
1415    ENDIF
1416
1417!
1418!-- Special actions during the run
1419    WRITE ( io, 470 )
1420    IF ( create_disturbances )  THEN
1421       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1422                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1423                          zu(disturbance_level_ind_t), disturbance_level_ind_t
1424       IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
1425          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1426       ELSE
1427          WRITE ( io, 473 )  disturbance_energy_limit
1428       ENDIF
1429       WRITE ( io, 474 )  TRIM( random_generator )
1430    ENDIF
1431    IF ( pt_surface_initial_change /= 0.0 )  THEN
1432       WRITE ( io, 475 )  pt_surface_initial_change
1433    ENDIF
[75]1434    IF ( humidity  .AND.  q_surface_initial_change /= 0.0 )  THEN
[1]1435       WRITE ( io, 476 )  q_surface_initial_change       
1436    ENDIF
1437    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0 )  THEN
1438       WRITE ( io, 477 )  q_surface_initial_change       
1439    ENDIF
1440
[60]1441    IF ( particle_advection )  THEN
[1]1442!
[60]1443!--    Particle attributes
1444       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1445                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
[117]1446                          end_time_prel, dt_sort_particles
[60]1447       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1448       IF ( random_start_position )  WRITE ( io, 481 )
1449       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1450       WRITE ( io, 495 )  total_number_of_particles
1451       IF ( maximum_number_of_tailpoints /= 0 )  THEN
1452          WRITE ( io, 483 )  maximum_number_of_tailpoints
1453          IF ( minimum_tailpoint_distance /= 0 )  THEN
1454             WRITE ( io, 484 )  total_number_of_tails,      &
1455                                minimum_tailpoint_distance, &
1456                                maximum_tailpoint_age
1457          ENDIF
[1]1458       ENDIF
[60]1459       IF ( dt_write_particle_data /= 9999999.9 )  THEN
1460          WRITE ( io, 485 )  dt_write_particle_data
1461          output_format = ''
1462          IF ( netcdf_output )  THEN
1463             IF ( netcdf_64bit )  THEN
1464                output_format = 'netcdf (64 bit offset) and binary'
1465             ELSE
1466                output_format = 'netcdf and binary'
1467             ENDIF
[1]1468          ELSE
[60]1469             output_format = 'binary'
[1]1470          ENDIF
[292]1471          WRITE ( io, 344 )  output_format
[1]1472       ENDIF
[60]1473       IF ( dt_dopts /= 9999999.9 )  WRITE ( io, 494 )  dt_dopts
1474       IF ( write_particle_statistics )  WRITE ( io, 486 )
[1]1475
[60]1476       WRITE ( io, 487 )  number_of_particle_groups
[1]1477
[60]1478       DO  i = 1, number_of_particle_groups
1479          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9 )  THEN
1480             WRITE ( io, 490 )  i, 0.0
1481             WRITE ( io, 492 )
[1]1482          ELSE
[60]1483             WRITE ( io, 490 )  i, radius(i)
1484             IF ( density_ratio(i) /= 0.0 )  THEN
1485                WRITE ( io, 491 )  density_ratio(i)
1486             ELSE
1487                WRITE ( io, 492 )
1488             ENDIF
[1]1489          ENDIF
[60]1490          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1491                             pdx(i), pdy(i), pdz(i)
[336]1492          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
[60]1493       ENDDO
[1]1494
[60]1495    ENDIF
[1]1496
[60]1497
[1]1498!
1499!-- Parameters of 1D-model
1500    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1501       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1502                          mixing_length_1d, dissipation_1d
1503       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1504          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1505       ENDIF
1506    ENDIF
1507
1508!
1509!-- User-defined informations
1510    CALL user_header( io )
1511
1512    WRITE ( io, 99 )
1513
1514!
1515!-- Write buffer contents to disc immediately
[82]1516    CALL local_flush( io )
[1]1517
1518!
1519!-- Here the FORMATs start
1520
1521 99 FORMAT (1X,78('-'))
[200]1522100 FORMAT (/1X,'***************************',9X,42('-')/        &
1523            1X,'* ',A,' *',9X,A/                               &
1524            1X,'***************************',9X,42('-'))
[291]1525101 FORMAT (37X,'coupled run using MPI-',I1,': ',A/ &
[102]1526            37X,42('-'))
[200]1527102 FORMAT (/' Date:              ',A8,9X,'Run:       ',A20/      &
1528            ' Time:              ',A8,9X,'Run-No.:   ',I2.2/     &
1529            ' Run on host:     ',A10)
[1]1530#if defined( __parallel )
[200]1531103 FORMAT (' Number of PEs:',8X,I5,9X,'Processor grid (x,y): (',I3,',',I3, &
[1]1532              ')',1X,A)
[200]1533104 FORMAT (' Number of PEs:',8X,I5,9X,'Tasks:',I4,'   threads per task:',I4/ &
[1]1534              37X,'Processor grid (x,y): (',I3,',',I3,')',1X,A)
[102]1535105 FORMAT (37X,'One additional PE is used to handle'/37X,'the dvrp output!')
1536106 FORMAT (37X,'A 1d-decomposition along x is forced'/ &
[1]1537            37X,'because the job is running on an SMP-cluster')
[102]1538107 FORMAT (37X,'A 1d-decomposition along ',A,' is used')
[1]1539#endif
1540110 FORMAT (/' Numerical Schemes:'/ &
1541             ' -----------------'/)
1542111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1543112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1544            '     Iterations (initial/other): ',I3,'/',I3,'  omega = ',F5.3)
1545113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1546                  ' or Upstream')
1547114 FORMAT (' --> Momentum advection via Upstream-Spline-Scheme')
1548115 FORMAT ('     Tendencies are smoothed via Long-Filter with factor ',F5.3) 
1549116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1550                  ' or Upstream')
1551117 FORMAT (' --> Scalar advection via Upstream-Spline-Scheme')
1552118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1553119 FORMAT (' --> Galilei-Transform applied to horizontal advection', &
1554            '     Translation velocity = ',A/ &
1555            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1556120 FORMAT (' --> Time differencing scheme: leapfrog only (no euler in case', &
1557                  ' of timestep changes)')
1558121 FORMAT (' --> Time differencing scheme: leapfrog + euler in case of', &
1559                  ' timestep changes')
1560122 FORMAT (' --> Time differencing scheme: ',A)
[108]1561123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
[1]1562            '     maximum damping coefficient: ',F5.3, ' 1/s')
1563124 FORMAT ('     Spline-overshoots are being suppressed')
1564125 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1565                  ' of'/                                                       &
1566            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1567126 FORMAT ('     Upstream-Scheme is used if Upstream-differences fall short', &
1568                  ' of'/                                                       &
1569            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1570127 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1571            '     delta_u = ',F6.4,4X,'delta_v = ',F6.4,4X,'delta_w = ',F6.4)
1572128 FORMAT ('     The following absolute overshoot differences are tolerated:'/&
1573            '     delta_e = ',F6.4,4X,'delta_pt = ',F6.4)
1574129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1575130 FORMAT (' --> Additional prognostic equation for the total water content')
1576131 FORMAT (' --> Parameterization of condensation processes via (0%-or100%)')
1577132 FORMAT (' --> Parameterization of long-wave radiation processes via'/ &
1578            '     effective emissivity scheme')
1579133 FORMAT (' --> Precipitation parameterization via Kessler-Scheme')
1580134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1581135 FORMAT (' --> Solve perturbation pressure via multigrid method (', &
1582                  A,'-cycle)'/ &
1583            '     number of grid levels:                   ',I2/ &
1584            '     Gauss-Seidel red/black iterations:       ',I2)
1585136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1586                  I3,')')
1587137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1588            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1589                  I3,')'/ &
1590            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1591                  I3,')')
1592138 FORMAT ('     Using hybrid version for 1d-domain-decomposition')
[63]1593139 FORMAT (' --> Loop optimization method: ',A)
[1]1594140 FORMAT ('     maximum residual allowed:                ',E10.3)
1595141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1596142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1597                  'step')
[87]1598143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1599                  'kinetic energy')
[1]1600150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
[241]1601                  'conserved'/ &
1602            '     using the ',A,' mode')
1603151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
[306]1604152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
1605           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
1606           /'     starting from dp_level_b =', F8.3, 'm', A /)
[411]1607153 FORMAT (' --> Large-scale vertical motion is used in the ', &
1608                  'prognostic equation for')
1609154 FORMAT ('     the potential temperature')
[1]1610200 FORMAT (//' Run time and time step information:'/ &
1611             ' ----------------------------------'/)
1612201 FORMAT ( ' Timestep:          variable     maximum value: ',F6.3,' s', &
1613             '    CFL-factor: ',F4.2)
1614202 FORMAT ( ' Timestep:       dt = ',F6.3,' s'/)
1615203 FORMAT ( ' Start time:       ',F9.3,' s'/ &
1616             ' End time:         ',F9.3,' s')
1617204 FORMAT ( A,F9.3,' s')
1618205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
1619206 FORMAT (/' Time reached:     ',F9.3,' s'/ &
1620             ' CPU-time used:    ',F9.3,' s     per timestep:               ', &
1621               '  ',F9.3,' s'/                                                 &
1622             '                                   per second of simulated tim', &
1623               'e: ',F9.3,' s')
[291]1624207 FORMAT ( A/' Coupling start time:',F9.3,' s')
[1]1625250 FORMAT (//' Computational grid and domain size:'/ &
1626              ' ----------------------------------'// &
1627              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1628              ' m    dz =    ',F7.3,' m'/ &
1629              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1630              ' m  z(u) = ',F10.3,' m'/)
1631252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
1632              ' factor: ',F5.3/ &
1633            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1634254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1635            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1636255 FORMAT (' Subdomains have equal size')
1637256 FORMAT (' Subdomains at the upper edges of the virtual processor grid ', &
1638              'have smaller sizes'/                                          &
1639            ' Size of smallest subdomain:    (  ',I4,',   ',I4,',   ',I4,')')
1640260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1641             ' degrees')
1642270 FORMAT (//' Topography informations:'/ &
1643              ' -----------------------'// &
1644              1X,'Topography: ',A)
1645271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1646              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1647                ' / ',I4)
[240]1648272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
1649              ' direction' / &
1650              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
1651              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
[256]1652278 FORMAT (' Topography grid definition convention:'/ &
1653            ' cell edge (staggered grid points'/  &
1654            ' (u in x-direction, v in y-direction))' /)
1655279 FORMAT (' Topography grid definition convention:'/ &
1656            ' cell center (scalar grid points)' /)
[138]1657280 FORMAT (//' Vegetation canopy (drag) model:'/ &
1658              ' ------------------------------'// &
1659              ' Canopy mode: ', A / &
1660              ' Canopy top: ',I4 / &
1661              ' Leaf drag coefficient: ',F6.2 /)
[153]1662281 FORMAT (/ ' Scalar_exchange_coefficient: ',F6.2 / &
1663              ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /)
1664282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s')
1665283 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
[138]1666              ' Height:              ',A,'  m'/ &
1667              ' Leaf area density:   ',A,'  m**2/m**3'/ &
1668              ' Gradient:            ',A,'  m**2/m**4'/ &
1669              ' Gridpoint:           ',A)
1670               
[1]1671300 FORMAT (//' Boundary conditions:'/ &
1672             ' -------------------'// &
1673             '                     p                    uv             ', &
1674             '                   pt'// &
1675             ' B. bound.: ',A/ &
1676             ' T. bound.: ',A)
[97]1677301 FORMAT (/'                     ',A// &
[1]1678             ' B. bound.: ',A/ &
1679             ' T. bound.: ',A)
[19]1680303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1681304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1682305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
1683               'computational u,v-level:'// &
[1]1684             '       zp = ',F6.2,' m   z0 = ',F6.4,' m   kappa = ',F4.2/ &
1685             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
[97]1686306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
[1]1687307 FORMAT ('       Heatflux has a random normal distribution')
1688308 FORMAT ('       Predefined surface temperature')
[97]1689309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
[1]1690310 FORMAT (//'    1D-Model:'// &
1691             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1692311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1693312 FORMAT ('       Predefined surface humidity')
1694313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1695314 FORMAT ('       Predefined scalar value at the surface')
[19]1696315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
[102]1697316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1698                    'atmosphere model')
[1]1699317 FORMAT (//' Lateral boundaries:'/ &
1700            '       left/right:  ',A/    &
1701            '       north/south: ',A)
1702318 FORMAT (/'       outflow damping layer width: ',I3,' gridpoints with km_', &
1703                    'max =',F5.1,' m**2/s')
[151]1704319 FORMAT ('       turbulence recycling at inflow switched on'/ &
1705            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1706            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
1707320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
[103]1708            '                                          v: ',F9.6,' m**2/s**2')
[151]1709325 FORMAT (//' List output:'/ &
[1]1710             ' -----------'//  &
1711            '    1D-Profiles:'/    &
1712            '       Output every             ',F8.2,' s')
[151]1713326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
[1]1714            '       Averaging input every    ',F8.2,' s')
1715330 FORMAT (//' Data output:'/ &
1716             ' -----------'/)
1717331 FORMAT (/'    1D-Profiles:')
1718332 FORMAT (/'       ',A)
1719333 FORMAT ('       Output every             ',F8.2,' s',/ &
1720            '       Time averaged over       ',F8.2,' s'/ &
1721            '       Averaging input every    ',F8.2,' s')
1722334 FORMAT (/'    2D-Arrays',A,':')
1723335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1724            '       Output every             ',F8.2,' s  ',A/ &
1725            '       Cross sections at ',A1,' = ',A/ &
1726            '       scalar-coordinates:   ',A,' m'/)
1727336 FORMAT (/'    3D-Arrays',A,':')
1728337 FORMAT (/'       Arrays: ',A/ &
1729            '       Output every             ',F8.2,' s  ',A/ &
1730            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1731338 FORMAT ('       Compressed data output'/ &
1732            '       Decimal precision: ',A/)
1733339 FORMAT ('       No output during initial ',F8.2,' s')
1734340 FORMAT (/'    Time series:')
1735341 FORMAT ('       Output every             ',F8.2,' s'/)
1736342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1737            '       Output every             ',F8.2,' s  ',A/ &
1738            '       Time averaged over       ',F8.2,' s'/ &
1739            '       Averaging input every    ',F8.2,' s'/ &
1740            '       Cross sections at ',A1,' = ',A/ &
1741            '       scalar-coordinates:   ',A,' m'/)
1742343 FORMAT (/'       Arrays: ',A/ &
1743            '       Output every             ',F8.2,' s  ',A/ &
1744            '       Time averaged over       ',F8.2,' s'/ &
1745            '       Averaging input every    ',F8.2,' s'/ &
1746            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
[292]1747344 FORMAT ('       Output format: ',A/)
[410]1748345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
1749            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
1750            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
1751            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
1752346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
1753347 FORMAT ('       Variables: ',A/ &
1754            '       Output every             ',F8.2,' s')
1755348 FORMAT ('       Variables: ',A/ &
1756            '       Output every             ',F8.2,' s'/ &
1757            '       Time averaged over       ',F8.2,' s'/ &
1758            '       Averaging input every    ',F8.2,' s')
1759349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1760            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
1761            13('       ',8(F8.2,',')/) )
1762350 FORMAT (/'       Output locations in ',A,'-direction: ', &
1763            'all gridpoints along ',A,'-direction (default).' )
1764351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1765            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
1766            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
[1]1767#if defined( __dvrp_graphics )
1768360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
1769            '       Output every      ',F7.1,' s'/ &
1770            '       Output mode:      ',A/ &
1771            '       Host / User:      ',A,' / ',A/ &
1772            '       Directory:        ',A// &
1773            '       The sequence contains:')
[337]1774361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
1775            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
1776362 FORMAT (/'       Slicer plane ',A/ &
[336]1777            '       Slicer limits: [',F6.2,',',F6.2,']')
[337]1778363 FORMAT (/'       Particles'/ &
[336]1779            '          particle size:  ',F7.2,' m')
1780364 FORMAT ('          particle ',A,' controlled by "',A,'" with interval [', &
1781                       F6.2,',',F6.2,']')
[337]1782365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
[336]1783            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
1784                     ')'/ &
1785            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
1786            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
[337]1787366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
[336]1788367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
[1]1789#endif
1790#if defined( __spectra )
1791370 FORMAT ('    Spectra:')
1792371 FORMAT ('       Output every ',F7.1,' s'/)
1793372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
1794            '       Directions: ', 10(A5,',')/                         &
[189]1795            '       height levels  k = ', 20(I3,',')/                  &
1796            '                          ', 20(I3,',')/                  &
1797            '                          ', 20(I3,',')/                  &
1798            '                          ', 20(I3,',')/                  &
1799            '                          ', 19(I3,','),I3,'.'/           &
[1]1800            '       height levels selected for standard plot:'/        &
[189]1801            '                      k = ', 20(I3,',')/                  &
1802            '                          ', 20(I3,',')/                  &
1803            '                          ', 20(I3,',')/                  &
1804            '                          ', 20(I3,',')/                  &
1805            '                          ', 19(I3,','),I3,'.'/           &
[1]1806            '       Time averaged over ', F7.1, ' s,' /                &
1807            '       Profiles for the time averaging are taken every ', &
1808                    F6.1,' s')
1809#endif
1810400 FORMAT (//' Physical quantities:'/ &
1811              ' -------------------'/)
1812410 FORMAT ('    Angular velocity    :   omega = ',E9.3,' rad/s'/  &
1813            '    Geograph. latitude  :   phi   = ',F4.1,' degr'/   &
1814            '    Coriolis parameter  :   f     = ',F9.6,' 1/s'/    &
1815            '                            f*    = ',F9.6,' 1/s')
1816411 FORMAT (/'    Gravity             :   g     = ',F4.1,' m/s**2')
[97]1817412 FORMAT (/'    Reference density in buoyancy terms: ',F8.3,' kg/m**3')
1818413 FORMAT (/'    Reference temperature in buoyancy terms: ',F8.4,' K')
[57]1819415 FORMAT (/'    Cloud physics parameters:'/ &
[1]1820             '    ------------------------'/)
[57]1821416 FORMAT ('        Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
[1]1822            '        Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
1823            '        Density of air     :   rho_0 = ',F5.3,' kg/m**3'/  &
1824            '        Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
1825            '        Vapourization heat :   L_v   = ',E8.2,' J/kg')
1826420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1827            '       Height:        ',A,'  m'/ &
1828            '       Temperature:   ',A,'  K'/ &
1829            '       Gradient:      ',A,'  K/100m'/ &
1830            '       Gridpoint:     ',A)
1831421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1832            '       Height:      ',A,'  m'/ &
1833            '       Humidity:    ',A,'  kg/kg'/ &
1834            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1835            '       Gridpoint:   ',A)
1836422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1837            '       Height:                  ',A,'  m'/ &
1838            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1839            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1840            '       Gridpoint:               ',A)
1841423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1842            '       Height:      ',A,'  m'/ &
1843            '       ug:          ',A,'  m/s'/ &
1844            '       Gradient:    ',A,'  1/100s'/ &
1845            '       Gridpoint:   ',A)
1846424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1847            '       Height:      ',A,'  m'/ &
[97]1848            '       vg:          ',A,'  m/s'/ &
[1]1849            '       Gradient:    ',A,'  1/100s'/ &
1850            '       Gridpoint:   ',A)
[97]1851425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
1852            '       Height:     ',A,'  m'/ &
1853            '       Salinity:   ',A,'  psu'/ &
1854            '       Gradient:   ',A,'  psu/100m'/ &
1855            '       Gridpoint:  ',A)
[411]1856426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
1857            '       Height:      ',A,'  m'/ &
1858            '       w_subs:      ',A,'  m/s'/ &
1859            '       Gradient:    ',A,'  (m/s)/100m'/ &
1860            '       Gridpoint:   ',A)
[1]1861450 FORMAT (//' LES / Turbulence quantities:'/ &
1862              ' ---------------------------'/)
1863451 FORMAT ('   Diffusion coefficients are constant:'/ &
1864            '   Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1865452 FORMAT ('   Mixing length is limited to the Prandtl mixing lenth.')
1866453 FORMAT ('   Mixing length is limited to ',F4.2,' * z')
1867454 FORMAT ('   TKE is not allowed to fall below ',E9.2,' (m/s)**2')
[108]1868455 FORMAT ('   initial TKE is prescribed as ',E9.2,' (m/s)**2')
[1]1869470 FORMAT (//' Actions during the simulation:'/ &
1870              ' -----------------------------'/)
[94]1871471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1872            '    Disturbance amplitude           :     ',F4.2, ' m/s'/       &
1873            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1874            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
[1]1875472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1876                 ' to i/j =',I4)
1877473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1878                 1X,F5.3, ' m**2/s**2')
1879474 FORMAT ('    Random number generator used    : ',A/)
1880475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1881                 'respectively, if'/ &
1882            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1883                 ' 3D-simulation'/)
1884476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1885                 'respectively, if the'/ &
1886            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1887                 ' the 3D-simulation'/)
1888477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1889                 'respectively, if the'/ &
1890            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1891                 ' the 3D-simulation'/)
1892480 FORMAT ('    Particles:'/ &
1893            '    ---------'// &
1894            '       Particle advection is active (switched on at t = ', F7.1, &
1895                    ' s)'/ &
1896            '       Start of new particle generations every  ',F6.1,' s'/ &
1897            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
1898            '                            bottom:     ', A, ' top:         ', A/&
1899            '       Maximum particle age:                 ',F9.1,' s'/ &
[117]1900            '       Advection stopped at t = ',F9.1,' s'/ &
1901            '       Particles are sorted every ',F9.1,' s'/)
[1]1902481 FORMAT ('       Particles have random start positions'/)
[336]1903482 FORMAT ('          Particles are advected only horizontally'/)
[1]1904483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
1905484 FORMAT ('            Number of tails of the total domain: ',I10/ &
1906            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
1907            '            Maximum age of the end of the tail:  ',F8.2,' s')
1908485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
1909486 FORMAT ('       Particle statistics are written on file'/)
1910487 FORMAT ('       Number of particle groups: ',I2/)
1911488 FORMAT ('       SGS velocity components are used for particle advection'/ &
1912            '          minimum timestep for advection: ', F7.5/)
1913489 FORMAT ('       Number of particles simultaneously released at each ', &
1914                    'point: ', I5/)
1915490 FORMAT ('       Particle group ',I2,':'/ &
1916            '          Particle radius: ',E10.3, 'm')
1917491 FORMAT ('          Particle inertia is activated'/ &
1918            '             density_ratio (rho_fluid/rho_particle) = ',F5.3/)
1919492 FORMAT ('          Particles are advected only passively (no inertia)'/)
1920493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
1921            '                                         y:',F8.1,' - ',F8.1,' m'/&
1922            '                                         z:',F8.1,' - ',F8.1,' m'/&
1923            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
1924                       ' m  dz = ',F8.1,' m'/)
1925494 FORMAT ('       Output of particle time series in NetCDF format every ', &
1926                    F8.2,' s'/)
1927495 FORMAT ('       Number of particles in total domain: ',I10/)
1928500 FORMAT (//' 1D-Model parameters:'/                           &
1929              ' -------------------'//                           &
1930            '    Simulation time:                   ',F8.1,' s'/ &
1931            '    Run-controll output every:         ',F8.1,' s'/ &
1932            '    Vertical profile output every:     ',F8.1,' s'/ &
1933            '    Mixing length calculation:         ',A/         &
1934            '    Dissipation calculation:           ',A/)
1935502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
1936
1937
1938 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.