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

Last change on this file since 1827 was 1827, checked in by maronga, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 87.3 KB
RevLine 
[1682]1!> @file header.f90
[1036]2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
[1818]16! Copyright 1997-2016 Leibniz Universitaet Hannover
[1036]17!--------------------------------------------------------------------------------!
18!
[254]19! Current revisions:
[1]20! -----------------
[1827]21!
22!
[1485]23! Former revisions:
24! -----------------
25! $Id: header.f90 1827 2016-04-07 12:12:23Z maronga $
26!
[1827]27! 1826 2016-04-07 12:01:39Z maronga
28! Moved radiation model header output to the respective module.
29! Moved canopy model header output to the respective module.
30!
[1823]31! 1822 2016-04-07 07:49:42Z hoffmann
32! Tails removed. icloud_scheme replaced by microphysics_*
33!
[1818]34! 1817 2016-04-06 15:44:20Z maronga
35! Moved land_surface_model header output to the respective module.
36!
[1809]37! 1808 2016-04-05 19:44:00Z raasch
38! routine local_flush replaced by FORTRAN statement
39!
[1798]40! 1797 2016-03-21 16:50:28Z raasch
41! output of nesting datatransfer mode
42!
[1792]43! 1791 2016-03-11 10:41:25Z raasch
44! output of nesting informations of all domains
45!
[1789]46! 1788 2016-03-10 11:01:04Z maronga
47! Parameter dewfall removed
48!
[1787]49! 1786 2016-03-08 05:49:27Z raasch
50! cpp-direktives for spectra removed
51!
[1784]52! 1783 2016-03-06 18:36:17Z raasch
53! netcdf module and variable names changed, output of netcdf_deflate
54!
[1765]55! 1764 2016-02-28 12:45:19Z raasch
56! output of nesting informations
57!
[1698]58! 1697 2015-10-28 17:14:10Z raasch
59! small E- and F-FORMAT changes to avoid informative compiler messages about
60! insufficient field width
61!
[1692]62! 1691 2015-10-26 16:17:44Z maronga
63! Renamed prandtl_layer to constant_flux_layer, renames rif_min/rif_max to
64! zeta_min/zeta_max.
65!
[1683]66! 1682 2015-10-07 23:56:08Z knoop
67! Code annotations made doxygen readable
68!
[1676]69! 1675 2015-10-02 08:28:59Z gronemeier
70! Bugfix: Definition of topography grid levels
71!
[1662]72! 1660 2015-09-21 08:15:16Z gronemeier
73! Bugfix: Definition of building/street canyon height if vertical grid stretching
74!         starts below the maximum topography height.
75!
[1591]76! 1590 2015-05-08 13:56:27Z maronga
77! Bugfix: Added TRIM statements for character strings for LSM and radiation code
78!
[1586]79! 1585 2015-04-30 07:05:52Z maronga
80! Further output for radiation model(s).
81!
[1576]82! 1575 2015-03-27 09:56:27Z raasch
83! adjustments for psolver-queries, output of seed_follows_topography
84!
[1561]85! 1560 2015-03-06 10:48:54Z keck
86! output for recycling y shift
87!
[1558]88! 1557 2015-03-05 16:43:04Z suehring
89! output for monotonic limiter
90!
[1552]91! 1551 2015-03-03 14:18:16Z maronga
92! Added informal output for land surface model and radiation model. Removed typo.
93!
[1497]94! 1496 2014-12-02 17:25:50Z maronga
95! Renamed: "radiation -> "cloud_top_radiation"
96!
[1485]97! 1484 2014-10-21 10:53:05Z kanani
[1484]98! Changes due to new module structure of the plant canopy model:
99!   module plant_canopy_model_mod and output for new canopy model parameters
100!   (alpha_lad, beta_lad, lai_beta,...) added,
101!   drag_coefficient, leaf_surface_concentration and scalar_exchange_coefficient
102!   renamed to canopy_drag_coeff, leaf_surface_conc and leaf_scalar_exch_coeff,
103!   learde renamed leaf_area_density.
104! Bugfix: DO-WHILE-loop for lad header information additionally restricted
105! by maximum number of gradient levels (currently 10)
[1483]106!
107! 1482 2014-10-18 12:34:45Z raasch
108! information about calculated or predefined virtual processor topology adjusted
109!
[1469]110! 1468 2014-09-24 14:06:57Z maronga
111! Adapted for use on up to 6-digit processor cores
112!
[1430]113! 1429 2014-07-15 12:53:45Z knoop
114! header exended to provide ensemble_member_nr if specified
115!
[1377]116! 1376 2014-04-26 11:21:22Z boeske
117! Correction of typos
118!
[1366]119! 1365 2014-04-22 15:03:56Z boeske
120! New section 'Large scale forcing and nudging':
121! output of large scale forcing and nudging information,
122! new section for initial profiles created
123!
[1360]124! 1359 2014-04-11 17:15:14Z hoffmann
125! dt_sort_particles removed
126!
[1354]127! 1353 2014-04-08 15:21:23Z heinze
128! REAL constants provided with KIND-attribute
129!
[1329]130! 1327 2014-03-21 11:00:16Z raasch
131! parts concerning iso2d and avs output removed,
132! -netcdf output queries
133!
[1325]134! 1324 2014-03-21 09:13:16Z suehring
135! Bugfix: module spectrum added
136!
[1323]137! 1322 2014-03-20 16:38:49Z raasch
138! REAL functions provided with KIND-attribute,
139! some REAL constants defined as wp-kind
140!
[1321]141! 1320 2014-03-20 08:40:49Z raasch
[1320]142! ONLY-attribute added to USE-statements,
143! kind-parameters added to all INTEGER and REAL declaration statements,
144! kinds are defined in new module kinds,
145! revision history before 2012 removed,
146! comment fields (!:) to be used for variable explanations added to
147! all variable declaration statements
[1321]148!
[1309]149! 1308 2014-03-13 14:58:42Z fricke
150! output of the fixed number of output time levels
151! output_format adjusted for masked data if netcdf_data_format > 5
152!
[1300]153! 1299 2014-03-06 13:15:21Z heinze
154! output for using large_scale subsidence in combination
155! with large_scale_forcing
156! reformatting, more detailed explanations
157!
[1242]158! 1241 2013-10-30 11:36:58Z heinze
159! output for nudging + large scale forcing from external file
160!
[1217]161! 1216 2013-08-26 09:31:42Z raasch
162! output for transpose_compute_overlap
163!
[1213]164! 1212 2013-08-15 08:46:27Z raasch
165! output for poisfft_hybrid removed
166!
[1182]167! 1179 2013-06-14 05:57:58Z raasch
168! output of reference_state, use_reference renamed use_single_reference_value
169!
[1160]170! 1159 2013-05-21 11:58:22Z fricke
171! +use_cmax
172!
[1116]173! 1115 2013-03-26 18:16:16Z hoffmann
174! descriptions for Seifert-Beheng-cloud-physics-scheme added
175!
[1112]176! 1111 2013-03-08 23:54:10Z raasch
177! output of accelerator board information
178! ibc_p_b = 2 removed
179!
[1109]180! 1108 2013-03-05 07:03:32Z raasch
181! bugfix for r1106
182!
[1107]183! 1106 2013-03-04 05:31:38Z raasch
184! some format changes for coupled runs
185!
[1093]186! 1092 2013-02-02 11:24:22Z raasch
187! unused variables removed
188!
[1037]189! 1036 2012-10-22 13:43:42Z raasch
190! code put under GPL (PALM 3.9)
191!
[1035]192! 1031 2012-10-19 14:35:30Z raasch
193! output of netCDF data format modified
194!
[1017]195! 1015 2012-09-27 09:23:24Z raasch
[1365]196! output of Adjustment of mixing length to the Prandtl mixing length at first
[1017]197! grid point above ground removed
198!
[1004]199! 1003 2012-09-14 14:35:53Z raasch
200! output of information about equal/unequal subdomain size removed
201!
[1002]202! 1001 2012-09-13 14:08:46Z raasch
203! all actions concerning leapfrog- and upstream-spline-scheme removed
204!
[979]205! 978 2012-08-09 08:28:32Z fricke
206! -km_damp_max, outflow_damping_width
207! +pt_damping_factor, pt_damping_width
208! +z0h
209!
[965]210! 964 2012-07-26 09:14:24Z raasch
211! output of profil-related quantities removed
212!
[941]213! 940 2012-07-09 14:31:00Z raasch
214! Output in case of simulations for pure neutral stratification (no pt-equation
215! solved)
216!
[928]217! 927 2012-06-06 19:15:04Z raasch
218! output of masking_method for mg-solver
219!
[869]220! 868 2012-03-28 12:21:07Z raasch
221! translation velocity in Galilean transformation changed to 0.6 * ug
222!
[834]223! 833 2012-02-22 08:55:55Z maronga
224! Adjusted format for leaf area density
225!
[829]226! 828 2012-02-21 12:00:36Z raasch
227! output of dissipation_classes + radius_classes
228!
[826]229! 825 2012-02-19 03:03:44Z raasch
230! Output of cloud physics parameters/quantities complemented and restructured
231!
[1]232! Revision 1.1  1997/08/11 06:17:20  raasch
233! Initial revision
234!
235!
236! Description:
237! ------------
[1764]238!> Writing a header with all important information about the current run.
[1682]239!> This subroutine is called three times, two times at the beginning
240!> (writing information on files RUN_CONTROL and HEADER) and one time at the
241!> end of the run, then writing additional information about CPU-usage on file
242!> header.
[411]243!-----------------------------------------------------------------------------!
[1682]244 SUBROUTINE header
245 
[1]246
[1320]247    USE arrays_3d,                                                             &
[1660]248        ONLY:  pt_init, qsws, q_init, sa_init, shf, ug, vg, w_subs, zu, zw
[1320]249       
[1]250    USE control_parameters
[1320]251       
252    USE cloud_parameters,                                                      &
253        ONLY:  cp, curvature_solution_effects, c_sedimentation,                &
254               limiter_sedimentation, l_v, nc_const, r_d, ventilation_effect
255       
256    USE cpulog,                                                                &
257        ONLY:  log_point_s
258       
259    USE dvrp_variables,                                                        &
260        ONLY:  use_seperate_pe_for_dvrp_output
261       
262    USE grid_variables,                                                        &
263        ONLY:  dx, dy
264       
265    USE indices,                                                               &
266        ONLY:  mg_loc_ind, nnx, nny, nnz, nx, ny, nxl_mg, nxr_mg, nyn_mg,      &
267               nys_mg, nzt, nzt_mg
268       
269    USE kinds
[1817]270 
[1551]271    USE land_surface_model_mod,                                                &
[1817]272        ONLY: land_surface, lsm_header
[1551]273 
[1320]274    USE model_1d,                                                              &
275        ONLY:  damp_level_ind_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
276       
[1783]277    USE netcdf_interface,                                                      &
278        ONLY:  netcdf_data_format, netcdf_data_format_string, netcdf_deflate
279
[1320]280    USE particle_attributes,                                                   &
281        ONLY:  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, collision_kernel,     &
282               density_ratio, dissipation_classes, dt_min_part, dt_prel,       &
[1359]283               dt_write_particle_data, end_time_prel,                          &
[1822]284               number_of_particle_groups, particle_advection,                  &
285               particle_advection_start,                                       &
[1320]286               particles_per_point, pdx, pdy, pdz,  psb, psl, psn, psr, pss,   &
287               pst, radius, radius_classes, random_start_position,             &
[1575]288               seed_follows_topography,                                        &
[1822]289               total_number_of_particles, use_sgs_for_particles,               &
[1320]290               vertical_particle_advection, write_particle_statistics
291       
[1]292    USE pegrid
[1484]293
294    USE plant_canopy_model_mod,                                                &
[1826]295        ONLY:  pcm_header, plant_canopy
[1551]296
[1791]297    USE pmc_handle_communicator,                                               &
298        ONLY:  pmc_get_model_info
299
[1764]300    USE pmc_interface,                                                         &
[1797]301        ONLY:  nested_run, nesting_datatransfer_mode, nesting_mode
[1764]302
[1551]303    USE radiation_model_mod,                                                   &
[1826]304        ONLY:  radiation, radiation_header
[1324]305   
306    USE spectrum,                                                              &
307        ONLY:  comp_spectra_level, data_output_sp, plot_spectra_level,         &
308               spectra_direction
[1]309
310    IMPLICIT NONE
311
[1682]312    CHARACTER (LEN=1)  ::  prec                !<
[1320]313   
[1682]314    CHARACTER (LEN=2)  ::  do2d_mode           !<
[1320]315   
[1682]316    CHARACTER (LEN=5)  ::  section_chr         !<
[1320]317   
[1682]318    CHARACTER (LEN=10) ::  coor_chr            !<
319    CHARACTER (LEN=10) ::  host_chr            !<
[1320]320   
[1682]321    CHARACTER (LEN=16) ::  begin_chr           !<
[1320]322   
[1682]323    CHARACTER (LEN=26) ::  ver_rev             !<
[1791]324
325    CHARACTER (LEN=32) ::  cpl_name            !<
[1320]326   
[1682]327    CHARACTER (LEN=40) ::  output_format       !<
[1320]328   
[1682]329    CHARACTER (LEN=70) ::  char1               !<
330    CHARACTER (LEN=70) ::  char2               !<
331    CHARACTER (LEN=70) ::  dopr_chr            !<
332    CHARACTER (LEN=70) ::  do2d_xy             !<
333    CHARACTER (LEN=70) ::  do2d_xz             !<
334    CHARACTER (LEN=70) ::  do2d_yz             !<
335    CHARACTER (LEN=70) ::  do3d_chr            !<
336    CHARACTER (LEN=70) ::  domask_chr          !<
337    CHARACTER (LEN=70) ::  run_classification  !<
[1320]338   
[1826]339    CHARACTER (LEN=85) ::  r_upper             !<
340    CHARACTER (LEN=85) ::  r_lower             !<
[1320]341   
[1682]342    CHARACTER (LEN=86) ::  coordinates         !<
343    CHARACTER (LEN=86) ::  gradients           !<
344    CHARACTER (LEN=86) ::  slices              !<
345    CHARACTER (LEN=86) ::  temperatures        !<
346    CHARACTER (LEN=86) ::  ugcomponent         !<
347    CHARACTER (LEN=86) ::  vgcomponent         !<
[1]348
[1682]349    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)  !<
[410]350
[1791]351    INTEGER(iwp) ::  av             !<
352    INTEGER(iwp) ::  bh             !<
353    INTEGER(iwp) ::  blx            !<
354    INTEGER(iwp) ::  bly            !<
355    INTEGER(iwp) ::  bxl            !<
356    INTEGER(iwp) ::  bxr            !<
357    INTEGER(iwp) ::  byn            !<
358    INTEGER(iwp) ::  bys            !<
359    INTEGER(iwp) ::  ch             !<
360    INTEGER(iwp) ::  count          !<
361    INTEGER(iwp) ::  cpl_parent_id  !<
362    INTEGER(iwp) ::  cwx            !<
363    INTEGER(iwp) ::  cwy            !<
364    INTEGER(iwp) ::  cxl            !<
365    INTEGER(iwp) ::  cxr            !<
366    INTEGER(iwp) ::  cyn            !<
367    INTEGER(iwp) ::  cys            !<
368    INTEGER(iwp) ::  dim            !<
369    INTEGER(iwp) ::  i              !<
370    INTEGER(iwp) ::  io             !<
371    INTEGER(iwp) ::  j              !<
372    INTEGER(iwp) ::  k              !<
373    INTEGER(iwp) ::  l              !<
374    INTEGER(iwp) ::  ll             !<
375    INTEGER(iwp) ::  mpi_type       !<
376    INTEGER(iwp) ::  my_cpl_id      !<
377    INTEGER(iwp) ::  n              !<
378    INTEGER(iwp) ::  ncpl           !<
379    INTEGER(iwp) ::  npe_total      !<
[1320]380   
[1826]381
[1682]382    REAL(wp) ::  cpuseconds_per_simulated_second  !<
[1791]383    REAL(wp) ::  lower_left_coord_x               !< x-coordinate of nest domain
384    REAL(wp) ::  lower_left_coord_y               !< y-coordinate of nest domain
[1]385
386!
387!-- Open the output file. At the end of the simulation, output is directed
388!-- to unit 19.
389    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
390         .NOT. simulated_time_at_begin /= simulated_time )  THEN
391       io = 15   !  header output on file RUN_CONTROL
392    ELSE
393       io = 19   !  header output on file HEADER
394    ENDIF
395    CALL check_open( io )
396
397!
398!-- At the end of the run, output file (HEADER) will be rewritten with
[1551]399!-- new information
[1]400    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
401
402!
403!-- Determine kind of model run
404    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
[1764]405       run_classification = 'restart run'
[328]406    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[1764]407       run_classification = 'run with cyclic fill of 3D - prerun data'
[147]408    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
[1764]409       run_classification = 'run without 1D - prerun'
[197]410    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
[1764]411       run_classification = 'run with 1D - prerun'
[197]412    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
[1764]413       run_classification = 'run initialized by user'
[1]414    ELSE
[254]415       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
416       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
[1]417    ENDIF
[1764]418    IF ( nested_run )  run_classification = 'nested ' // run_classification
[97]419    IF ( ocean )  THEN
420       run_classification = 'ocean - ' // run_classification
421    ELSE
422       run_classification = 'atmosphere - ' // run_classification
423    ENDIF
[1]424
425!
426!-- Run-identification, date, time, host
427    host_chr = host(1:10)
[75]428    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
[102]429    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
[291]430    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
431#if defined( __mpi2 )
432       mpi_type = 2
433#else
434       mpi_type = 1
435#endif
436       WRITE ( io, 101 )  mpi_type, coupling_mode
437    ENDIF
[1108]438#if defined( __parallel )
[1353]439    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]440       IF ( coupling_start_time > simulated_time_at_begin )  THEN
441          WRITE ( io, 109 )
442       ELSE
443          WRITE ( io, 114 )
444       ENDIF
445    ENDIF
[1108]446#endif
[1429]447    IF ( ensemble_member_nr /= 0 )  THEN
448       WRITE ( io, 512 )  run_date, run_identifier, run_time, runnr,           &
449                       ADJUSTR( host_chr ), ensemble_member_nr
450    ELSE
451       WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr,           &
[102]452                       ADJUSTR( host_chr )
[1429]453    ENDIF
[1]454#if defined( __parallel )
[1482]455    IF ( npex == -1  .AND.  npey == -1 )  THEN
[1]456       char1 = 'calculated'
457    ELSE
458       char1 = 'predefined'
459    ENDIF
460    IF ( threads_per_task == 1 )  THEN
[102]461       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
[1]462    ELSE
[102]463       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
[1]464                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
465    ENDIF
[1111]466    IF ( num_acc_per_node /= 0 )  WRITE ( io, 117 )  num_acc_per_node   
[1]467    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
468           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
469         npex == -1  .AND.  pdims(2) == 1 )                      &
470    THEN
[102]471       WRITE ( io, 106 )
[1]472    ELSEIF ( pdims(2) == 1 )  THEN
[102]473       WRITE ( io, 107 )  'x'
[1]474    ELSEIF ( pdims(1) == 1 )  THEN
[102]475       WRITE ( io, 107 )  'y'
[1]476    ENDIF
[102]477    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
[759]478    IF ( numprocs /= maximum_parallel_io_streams )  THEN
479       WRITE ( io, 108 )  maximum_parallel_io_streams
480    ENDIF
[1111]481#else
482    IF ( num_acc_per_node /= 0 )  WRITE ( io, 120 )  num_acc_per_node
[1]483#endif
[1764]484
485!
486!-- Nesting informations
487    IF ( nested_run )  THEN
[1791]488
[1797]489       WRITE ( io, 600 )  TRIM( nesting_mode ),                                &
490                          TRIM( nesting_datatransfer_mode )
[1791]491       CALL pmc_get_model_info( ncpl = ncpl, cpl_id = my_cpl_id )
492
493       DO  n = 1, ncpl
494          CALL pmc_get_model_info( request_for_cpl_id = n, cpl_name = cpl_name,&
495                                   cpl_parent_id = cpl_parent_id,              &
496                                   lower_left_x = lower_left_coord_x,          &
497                                   lower_left_y = lower_left_coord_y,          &
498                                   npe_total = npe_total )
499          IF ( n == my_cpl_id )  THEN
500             char1 = '*'
501          ELSE
502             char1 = ' '
503          ENDIF
504          WRITE ( io, 601 )  TRIM( char1 ), n, cpl_parent_id, npe_total,       &
505                             lower_left_coord_x, lower_left_coord_y,           &
506                             TRIM( cpl_name )
507       ENDDO
[1764]508    ENDIF
[1]509    WRITE ( io, 99 )
510
511!
512!-- Numerical schemes
513    WRITE ( io, 110 )
514    IF ( psolver(1:7) == 'poisfft' )  THEN
515       WRITE ( io, 111 )  TRIM( fft_method )
[1216]516       IF ( transpose_compute_overlap )  WRITE( io, 115 )
[1]517    ELSEIF ( psolver == 'sor' )  THEN
518       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
[1575]519    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
520       WRITE ( io, 135 )  TRIM(psolver), cycle_mg, maximum_grid_level, ngsrb
[1]521       IF ( mg_cycles == -1 )  THEN
522          WRITE ( io, 140 )  residual_limit
523       ELSE
524          WRITE ( io, 141 )  mg_cycles
525       ENDIF
526       IF ( mg_switch_to_pe0_level == 0 )  THEN
527          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
528                             nzt_mg(1)
[197]529       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
[1]530          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
531                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
532                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
533                             nzt_mg(mg_switch_to_pe0_level),    &
534                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
535                             nzt_mg(1)
536       ENDIF
[927]537       IF ( masking_method )  WRITE ( io, 144 )
[1]538    ENDIF
539    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
540    THEN
541       WRITE ( io, 142 )
542    ENDIF
543
544    IF ( momentum_advec == 'pw-scheme' )  THEN
545       WRITE ( io, 113 )
[1299]546    ELSEIF (momentum_advec == 'ws-scheme' )  THEN
[667]547       WRITE ( io, 503 )
[1]548    ENDIF
549    IF ( scalar_advec == 'pw-scheme' )  THEN
550       WRITE ( io, 116 )
[667]551    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
552       WRITE ( io, 504 )
[1557]553    ELSEIF ( scalar_advec == 'ws-scheme-mono' )  THEN
554       WRITE ( io, 513 )
[1]555    ELSE
556       WRITE ( io, 118 )
557    ENDIF
[63]558
559    WRITE ( io, 139 )  TRIM( loop_optimization )
560
[1]561    IF ( galilei_transformation )  THEN
562       IF ( use_ug_for_galilei_tr )  THEN
[868]563          char1 = '0.6 * geostrophic wind'
[1]564       ELSE
565          char1 = 'mean wind in model domain'
566       ENDIF
567       IF ( simulated_time_at_begin == simulated_time )  THEN
568          char2 = 'at the start of the run'
569       ELSE
570          char2 = 'at the end of the run'
571       ENDIF
[1353]572       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ),                        &
573                          advected_distance_x/1000.0_wp,                       &
574                          advected_distance_y/1000.0_wp
[1]575    ENDIF
[1001]576    WRITE ( io, 122 )  timestep_scheme
[87]577    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
[1353]578    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
[108]579       IF ( .NOT. ocean )  THEN
580          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
581               rayleigh_damping_factor
582       ELSE
583          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
584               rayleigh_damping_factor
585       ENDIF
[1]586    ENDIF
[940]587    IF ( neutral )  WRITE ( io, 131 )  pt_surface
[75]588    IF ( humidity )  THEN
[1]589       IF ( .NOT. cloud_physics )  THEN
590          WRITE ( io, 129 )
591       ELSE
592          WRITE ( io, 130 )
593       ENDIF
594    ENDIF
595    IF ( passive_scalar )  WRITE ( io, 134 )
[240]596    IF ( conserve_volume_flow )  THEN
[241]597       WRITE ( io, 150 )  conserve_volume_flow_mode
598       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
599          WRITE ( io, 151 )  u_bulk, v_bulk
600       ENDIF
[240]601    ELSEIF ( dp_external )  THEN
602       IF ( dp_smooth )  THEN
[241]603          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
[240]604       ELSE
[241]605          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
[240]606       ENDIF
607    ENDIF
[1]608    WRITE ( io, 99 )
609
610!
[1551]611!-- Runtime and timestep information
[1]612    WRITE ( io, 200 )
613    IF ( .NOT. dt_fixed )  THEN
614       WRITE ( io, 201 )  dt_max, cfl_factor
615    ELSE
616       WRITE ( io, 202 )  dt
617    ENDIF
618    WRITE ( io, 203 )  simulated_time_at_begin, end_time
619
[1322]620    IF ( time_restart /= 9999999.9_wp  .AND. &
[1]621         simulated_time_at_begin == simulated_time )  THEN
[1322]622       IF ( dt_restart == 9999999.9_wp )  THEN
[1]623          WRITE ( io, 204 )  ' Restart at:       ',time_restart
624       ELSE
625          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
626       ENDIF
627    ENDIF
628
629    IF ( simulated_time_at_begin /= simulated_time )  THEN
630       i = MAX ( log_point_s(10)%counts, 1 )
[1353]631       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0_wp )  THEN
632          cpuseconds_per_simulated_second = 0.0_wp
[1]633       ELSE
634          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
635                                            ( simulated_time -    &
636                                              simulated_time_at_begin )
637       ENDIF
[1322]638       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum,      &
639                          log_point_s(10)%sum / REAL( i, KIND=wp ), &
[1]640                          cpuseconds_per_simulated_second
[1322]641       IF ( time_restart /= 9999999.9_wp  .AND.  time_restart < end_time )  THEN
642          IF ( dt_restart == 9999999.9_wp )  THEN
[1106]643             WRITE ( io, 204 )  ' Next restart at:     ',time_restart
[1]644          ELSE
[1106]645             WRITE ( io, 205 )  ' Next restart at:     ',time_restart, dt_restart
[1]646          ENDIF
647       ENDIF
648    ENDIF
649
[1324]650
[1]651!
[291]652!-- Start time for coupled runs, if independent precursor runs for atmosphere
[1106]653!-- and ocean are used or have been used. In this case, coupling_start_time
654!-- defines the time when the coupling is switched on.
[1353]655    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]656       WRITE ( io, 207 )  coupling_start_time
[291]657    ENDIF
658
659!
[1]660!-- Computational grid
[94]661    IF ( .NOT. ocean )  THEN
662       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
663       IF ( dz_stretch_level_index < nzt+1 )  THEN
664          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
665                             dz_stretch_factor, dz_max
666       ENDIF
667    ELSE
668       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
669       IF ( dz_stretch_level_index > 0 )  THEN
670          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
671                             dz_stretch_factor, dz_max
672       ENDIF
[1]673    ENDIF
674    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
675                       MIN( nnz+2, nzt+2 )
676    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
677
678!
[1365]679!-- Large scale forcing and nudging
680    WRITE ( io, 160 )
681    IF ( large_scale_forcing )  THEN
682       WRITE ( io, 162 )
683       WRITE ( io, 163 )
684
685       IF ( large_scale_subsidence )  THEN
686          IF ( .NOT. use_subsidence_tendencies )  THEN
687             WRITE ( io, 164 )
688          ELSE
689             WRITE ( io, 165 )
690          ENDIF
691       ENDIF
692
693       IF ( bc_pt_b == 'dirichlet' )  THEN
694          WRITE ( io, 180 )
695       ELSEIF ( bc_pt_b == 'neumann' )  THEN
696          WRITE ( io, 181 )
697       ENDIF
698
699       IF ( bc_q_b == 'dirichlet' )  THEN
700          WRITE ( io, 182 )
701       ELSEIF ( bc_q_b == 'neumann' )  THEN
702          WRITE ( io, 183 )
703       ENDIF
704
705       WRITE ( io, 167 )
706       IF ( nudging )  THEN
707          WRITE ( io, 170 )
708       ENDIF
709    ELSE
710       WRITE ( io, 161 )
711       WRITE ( io, 171 )
712    ENDIF
713    IF ( large_scale_subsidence )  THEN
714       WRITE ( io, 168 )
715       WRITE ( io, 169 )
716    ENDIF
717
718!
719!-- Profile for the large scale vertial velocity
720!-- Building output strings, starting with surface value
721    IF ( large_scale_subsidence )  THEN
722       temperatures = '   0.0'
723       gradients = '------'
724       slices = '     0'
725       coordinates = '   0.0'
726       i = 1
727       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
728
729          WRITE (coor_chr,'(E10.2,7X)')  &
730                                w_subs(subs_vertical_gradient_level_i(i))
731          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
732
733          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
734          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
735
736          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
737          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
738
739          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
740          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
741
742          IF ( i == 10 )  THEN
743             EXIT
744          ELSE
745             i = i + 1
746          ENDIF
747
748       ENDDO
749
750 
751       IF ( .NOT. large_scale_forcing )  THEN
752          WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
753                             TRIM( gradients ), TRIM( slices )
754       ENDIF
755
756
757    ENDIF
758
759!-- Profile of the geostrophic wind (component ug)
760!-- Building output strings
761    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
762    gradients = '------'
763    slices = '     0'
764    coordinates = '   0.0'
765    i = 1
766    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
767     
768       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
769       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
770
771       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
772       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
773
774       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
775       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
776
777       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
778       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
779
780       IF ( i == 10 )  THEN
781          EXIT
782       ELSE
783          i = i + 1
784       ENDIF
785
786    ENDDO
787
788    IF ( .NOT. large_scale_forcing )  THEN
789       WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
790                          TRIM( gradients ), TRIM( slices )
791    ENDIF
792
793!-- Profile of the geostrophic wind (component vg)
794!-- Building output strings
795    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
796    gradients = '------'
797    slices = '     0'
798    coordinates = '   0.0'
799    i = 1
800    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
801
802       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
803       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
804
805       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
806       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
807
808       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
809       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
810
811       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
812       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
813
814       IF ( i == 10 )  THEN
815          EXIT
816       ELSE
817          i = i + 1
818       ENDIF
819 
820    ENDDO
821
822    IF ( .NOT. large_scale_forcing )  THEN
823       WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
824                          TRIM( gradients ), TRIM( slices )
825    ENDIF
826
827!
[1]828!-- Topography
829    WRITE ( io, 270 )  topography
830    SELECT CASE ( TRIM( topography ) )
831
832       CASE ( 'flat' )
833          ! no actions necessary
834
835       CASE ( 'single_building' )
836          blx = INT( building_length_x / dx )
837          bly = INT( building_length_y / dy )
[1675]838          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
839          IF ( ABS( zw(bh  ) - building_height ) == &
840               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
[1]841
[1322]842          IF ( building_wall_left == 9999999.9_wp )  THEN
[1]843             building_wall_left = ( nx + 1 - blx ) / 2 * dx
844          ENDIF
[1353]845          bxl = INT ( building_wall_left / dx + 0.5_wp )
[1]846          bxr = bxl + blx
847
[1322]848          IF ( building_wall_south == 9999999.9_wp )  THEN
[1]849             building_wall_south = ( ny + 1 - bly ) / 2 * dy
850          ENDIF
[1353]851          bys = INT ( building_wall_south / dy + 0.5_wp )
[1]852          byn = bys + bly
853
854          WRITE ( io, 271 )  building_length_x, building_length_y, &
855                             building_height, bxl, bxr, bys, byn
856
[240]857       CASE ( 'single_street_canyon' )
[1675]858          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
859          IF ( ABS( zw(ch  ) - canyon_height ) == &
860               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
[1322]861          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[240]862!
863!--          Street canyon in y direction
864             cwx = NINT( canyon_width_x / dx )
[1322]865             IF ( canyon_wall_left == 9999999.9_wp )  THEN
[240]866                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
867             ENDIF
868             cxl = NINT( canyon_wall_left / dx )
869             cxr = cxl + cwx
870             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
871
[1322]872          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[240]873!
874!--          Street canyon in x direction
875             cwy = NINT( canyon_width_y / dy )
[1322]876             IF ( canyon_wall_south == 9999999.9_wp )  THEN
[240]877                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
878             ENDIF
879             cys = NINT( canyon_wall_south / dy )
880             cyn = cys + cwy
881             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
882          ENDIF
883
[1]884    END SELECT
885
[256]886    IF ( TRIM( topography ) /= 'flat' )  THEN
887       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
888          IF ( TRIM( topography ) == 'single_building' .OR.  &
889               TRIM( topography ) == 'single_street_canyon' )  THEN
890             WRITE ( io, 278 )
891          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
892             WRITE ( io, 279 )
893          ENDIF
894       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
895          WRITE ( io, 278 )
896       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
897          WRITE ( io, 279 )
898       ENDIF
899    ENDIF
900
[1826]901    IF ( plant_canopy )  CALL pcm_header ( io )
[138]902
[1817]903    IF ( land_surface )  CALL lsm_header ( io )
[1484]904
[1826]905    IF ( radiation )  CALL radiation_header ( io )
[1551]906
907!
[1]908!-- Boundary conditions
909    IF ( ibc_p_b == 0 )  THEN
[1826]910       r_lower = 'p(0)     = 0      |'
[1]911    ELSEIF ( ibc_p_b == 1 )  THEN
[1826]912       r_lower = 'p(0)     = p(1)   |'
[1]913    ENDIF
914    IF ( ibc_p_t == 0 )  THEN
[1826]915       r_upper  = 'p(nzt+1) = 0      |'
[1]916    ELSE
[1826]917       r_upper  = 'p(nzt+1) = p(nzt) |'
[1]918    ENDIF
919
920    IF ( ibc_uv_b == 0 )  THEN
[1826]921       r_lower = TRIM( r_lower ) // ' uv(0)     = -uv(1)                |'
[1]922    ELSE
[1826]923       r_lower = TRIM( r_lower ) // ' uv(0)     = uv(1)                 |'
[1]924    ENDIF
[132]925    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
[1826]926       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = 0                     |'
[132]927    ELSEIF ( ibc_uv_t == 0 )  THEN
[1826]928       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
[1]929    ELSE
[1826]930       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = uv(nzt)               |'
[1]931    ENDIF
932
933    IF ( ibc_pt_b == 0 )  THEN
[1551]934       IF ( land_surface )  THEN
[1826]935          r_lower = TRIM( r_lower ) // ' pt(0)     = from soil model'
[1551]936       ELSE
[1826]937          r_lower = TRIM( r_lower ) // ' pt(0)     = pt_surface'
[1551]938       ENDIF
[102]939    ELSEIF ( ibc_pt_b == 1 )  THEN
[1826]940       r_lower = TRIM( r_lower ) // ' pt(0)     = pt(1)'
[102]941    ELSEIF ( ibc_pt_b == 2 )  THEN
[1826]942       r_lower = TRIM( r_lower ) // ' pt(0)     = from coupled model'
[1]943    ENDIF
944    IF ( ibc_pt_t == 0 )  THEN
[1826]945       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt_top'
[19]946    ELSEIF( ibc_pt_t == 1 )  THEN
[1826]947       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt)'
[19]948    ELSEIF( ibc_pt_t == 2 )  THEN
[1826]949       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
[667]950
[1]951    ENDIF
952
[1826]953    WRITE ( io, 300 )  r_lower, r_upper
[1]954
955    IF ( .NOT. constant_diffusion )  THEN
956       IF ( ibc_e_b == 1 )  THEN
[1826]957          r_lower = 'e(0)     = e(1)'
[1]958       ELSE
[1826]959          r_lower = 'e(0)     = e(1) = (u*/0.1)**2'
[1]960       ENDIF
[1826]961       r_upper = 'e(nzt+1) = e(nzt) = e(nzt-1)'
[1]962
[1826]963       WRITE ( io, 301 )  'e', r_lower, r_upper       
[1]964
965    ENDIF
966
[97]967    IF ( ocean )  THEN
[1826]968       r_lower = 'sa(0)    = sa(1)'
[97]969       IF ( ibc_sa_t == 0 )  THEN
[1826]970          r_upper =  'sa(nzt+1) = sa_surface'
[1]971       ELSE
[1826]972          r_upper =  'sa(nzt+1) = sa(nzt)'
[1]973       ENDIF
[1826]974       WRITE ( io, 301 ) 'sa', r_lower, r_upper
[97]975    ENDIF
[1]976
[97]977    IF ( humidity )  THEN
978       IF ( ibc_q_b == 0 )  THEN
[1551]979          IF ( land_surface )  THEN
[1826]980             r_lower = 'q(0)     = from soil model'
[1551]981          ELSE
[1826]982             r_lower = 'q(0)     = q_surface'
[1551]983          ENDIF
984
[97]985       ELSE
[1826]986          r_lower = 'q(0)     = q(1)'
[97]987       ENDIF
988       IF ( ibc_q_t == 0 )  THEN
[1826]989          r_upper =  'q(nzt)   = q_top'
[97]990       ELSE
[1826]991          r_upper =  'q(nzt)   = q(nzt-1) + dq/dz'
[97]992       ENDIF
[1826]993       WRITE ( io, 301 ) 'q', r_lower, r_upper
[97]994    ENDIF
[1]995
[97]996    IF ( passive_scalar )  THEN
997       IF ( ibc_q_b == 0 )  THEN
[1826]998          r_lower = 's(0)     = s_surface'
[97]999       ELSE
[1826]1000          r_lower = 's(0)     = s(1)'
[97]1001       ENDIF
1002       IF ( ibc_q_t == 0 )  THEN
[1826]1003          r_upper =  's(nzt)   = s_top'
[97]1004       ELSE
[1826]1005          r_upper =  's(nzt)   = s(nzt-1) + ds/dz'
[97]1006       ENDIF
[1826]1007       WRITE ( io, 301 ) 's', r_lower, r_upper
[1]1008    ENDIF
1009
1010    IF ( use_surface_fluxes )  THEN
1011       WRITE ( io, 303 )
1012       IF ( constant_heatflux )  THEN
[1299]1013          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[1241]1014             WRITE ( io, 306 )  shf(0,0)
1015          ELSE
1016             WRITE ( io, 306 )  surface_heatflux
1017          ENDIF
[1]1018          IF ( random_heatflux )  WRITE ( io, 307 )
1019       ENDIF
[75]1020       IF ( humidity  .AND.  constant_waterflux )  THEN
[1299]1021          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[1241]1022             WRITE ( io, 311 ) qsws(0,0)
1023          ELSE
1024             WRITE ( io, 311 ) surface_waterflux
1025          ENDIF
[1]1026       ENDIF
1027       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
1028          WRITE ( io, 313 ) surface_waterflux
1029       ENDIF
1030    ENDIF
1031
[19]1032    IF ( use_top_fluxes )  THEN
1033       WRITE ( io, 304 )
[102]1034       IF ( coupling_mode == 'uncoupled' )  THEN
[151]1035          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
[102]1036          IF ( constant_top_heatflux )  THEN
1037             WRITE ( io, 306 )  top_heatflux
1038          ENDIF
1039       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1040          WRITE ( io, 316 )
[19]1041       ENDIF
[97]1042       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
1043          WRITE ( io, 309 )  top_salinityflux
1044       ENDIF
[75]1045       IF ( humidity  .OR.  passive_scalar )  THEN
[19]1046          WRITE ( io, 315 )
1047       ENDIF
1048    ENDIF
1049
[1691]1050    IF ( constant_flux_layer )  THEN
1051       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length,                     &
1052                          z0h_factor*roughness_length, kappa,                  &
1053                          zeta_min, zeta_max
[1]1054       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
[75]1055       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
[1]1056          WRITE ( io, 312 )
1057       ENDIF
1058       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
1059          WRITE ( io, 314 )
1060       ENDIF
1061    ELSE
1062       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
[1691]1063          WRITE ( io, 310 )  zeta_min, zeta_max
[1]1064       ENDIF
1065    ENDIF
1066
1067    WRITE ( io, 317 )  bc_lr, bc_ns
[707]1068    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1159]1069       WRITE ( io, 318 )  use_cmax, pt_damping_width, pt_damping_factor       
[151]1070       IF ( turbulent_inflow )  THEN
[1560]1071          IF ( .NOT. recycling_yshift ) THEN
1072             WRITE ( io, 319 )  recycling_width, recycling_plane, &
1073                                inflow_damping_height, inflow_damping_width
1074          ELSE
1075             WRITE ( io, 322 )  recycling_width, recycling_plane, &
1076                                inflow_damping_height, inflow_damping_width
1077          END IF
[151]1078       ENDIF
[1]1079    ENDIF
1080
1081!
[1365]1082!-- Initial Profiles
1083    WRITE ( io, 321 )
1084!
1085!-- Initial wind profiles
1086    IF ( u_profile(1) /= 9999999.9_wp )  WRITE ( io, 427 )
1087
1088!
1089!-- Initial temperature profile
1090!-- Building output strings, starting with surface temperature
1091    WRITE ( temperatures, '(F6.2)' )  pt_surface
1092    gradients = '------'
1093    slices = '     0'
1094    coordinates = '   0.0'
1095    i = 1
1096    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1097
1098       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1099       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1100
1101       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1102       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1103
1104       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1105       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1106
1107       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1108       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1109
1110       IF ( i == 10 )  THEN
1111          EXIT
1112       ELSE
1113          i = i + 1
1114       ENDIF
1115
1116    ENDDO
1117
1118    IF ( .NOT. nudging )  THEN
1119       WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1120                          TRIM( gradients ), TRIM( slices )
1121    ELSE
1122       WRITE ( io, 428 ) 
1123    ENDIF
1124
1125!
1126!-- Initial humidity profile
1127!-- Building output strings, starting with surface humidity
1128    IF ( humidity  .OR.  passive_scalar )  THEN
1129       WRITE ( temperatures, '(E8.1)' )  q_surface
1130       gradients = '--------'
1131       slices = '       0'
1132       coordinates = '     0.0'
1133       i = 1
1134       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1135         
1136          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1137          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1138
1139          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1140          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1141         
1142          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1143          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1144         
1145          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1146          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1147
1148          IF ( i == 10 )  THEN
1149             EXIT
1150          ELSE
1151             i = i + 1
1152          ENDIF
1153
1154       ENDDO
1155
1156       IF ( humidity )  THEN
1157          IF ( .NOT. nudging )  THEN
1158             WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1159                                TRIM( gradients ), TRIM( slices )
1160          ENDIF
1161       ELSE
1162          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1163                             TRIM( gradients ), TRIM( slices )
1164       ENDIF
1165    ENDIF
1166
1167!
1168!-- Initial salinity profile
1169!-- Building output strings, starting with surface salinity
1170    IF ( ocean )  THEN
1171       WRITE ( temperatures, '(F6.2)' )  sa_surface
1172       gradients = '------'
1173       slices = '     0'
1174       coordinates = '   0.0'
1175       i = 1
1176       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1177
1178          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1179          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1180
1181          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1182          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1183
1184          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1185          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1186
1187          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1188          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1189
1190          IF ( i == 10 )  THEN
1191             EXIT
1192          ELSE
1193             i = i + 1
1194          ENDIF
1195
1196       ENDDO
1197
1198       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1199                          TRIM( gradients ), TRIM( slices )
1200    ENDIF
1201
1202
1203!
[1]1204!-- Listing of 1D-profiles
[151]1205    WRITE ( io, 325 )  dt_dopr_listing
[1353]1206    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1207       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1208    ENDIF
1209
1210!
1211!-- DATA output
1212    WRITE ( io, 330 )
[1353]1213    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1214       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1215    ENDIF
1216
1217!
1218!-- 1D-profiles
[346]1219    dopr_chr = 'Profile:'
[1]1220    IF ( dopr_n /= 0 )  THEN
1221       WRITE ( io, 331 )
1222
1223       output_format = ''
[1783]1224       output_format = netcdf_data_format_string
1225       IF ( netcdf_deflate == 0 )  THEN
1226          WRITE ( io, 344 )  output_format
1227       ELSE
1228          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1229       ENDIF
[1]1230
1231       DO  i = 1, dopr_n
1232          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
1233          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
1234             WRITE ( io, 332 )  dopr_chr
1235             dopr_chr = '       :'
1236          ENDIF
1237       ENDDO
1238
1239       IF ( dopr_chr /= '' )  THEN
1240          WRITE ( io, 332 )  dopr_chr
1241       ENDIF
1242       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
[1353]1243       IF ( skip_time_dopr /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dopr
[1]1244    ENDIF
1245
1246!
1247!-- 2D-arrays
1248    DO  av = 0, 1
1249
1250       i = 1
1251       do2d_xy = ''
1252       do2d_xz = ''
1253       do2d_yz = ''
1254       DO  WHILE ( do2d(av,i) /= ' ' )
1255
1256          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
1257          do2d_mode = do2d(av,i)(l-1:l)
1258
1259          SELECT CASE ( do2d_mode )
1260             CASE ( 'xy' )
1261                ll = LEN_TRIM( do2d_xy )
1262                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1263             CASE ( 'xz' )
1264                ll = LEN_TRIM( do2d_xz )
1265                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1266             CASE ( 'yz' )
1267                ll = LEN_TRIM( do2d_yz )
1268                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1269          END SELECT
1270
1271          i = i + 1
1272
1273       ENDDO
1274
1275       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
1276              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
[1327]1277              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) ) )  THEN
[1]1278
1279          IF (  av == 0 )  THEN
1280             WRITE ( io, 334 )  ''
1281          ELSE
1282             WRITE ( io, 334 )  '(time-averaged)'
1283          ENDIF
1284
1285          IF ( do2d_at_begin )  THEN
1286             begin_chr = 'and at the start'
1287          ELSE
1288             begin_chr = ''
1289          ENDIF
1290
1291          output_format = ''
[1783]1292          output_format = netcdf_data_format_string
1293          IF ( netcdf_deflate == 0 )  THEN
1294             WRITE ( io, 344 )  output_format
1295          ELSE
1296             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1297          ENDIF
[1]1298
1299          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
1300             i = 1
1301             slices = '/'
1302             coordinates = '/'
1303!
[1551]1304!--          Building strings with index and coordinate information of the
[1]1305!--          slices
1306             DO  WHILE ( section(i,1) /= -9999 )
1307
1308                WRITE (section_chr,'(I5)')  section(i,1)
1309                section_chr = ADJUSTL( section_chr )
1310                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1311
[206]1312                IF ( section(i,1) == -1 )  THEN
[1353]1313                   WRITE (coor_chr,'(F10.1)')  -1.0_wp
[206]1314                ELSE
1315                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
1316                ENDIF
[1]1317                coor_chr = ADJUSTL( coor_chr )
1318                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1319
1320                i = i + 1
1321             ENDDO
1322             IF ( av == 0 )  THEN
1323                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
1324                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
1325                                   TRIM( coordinates )
[1353]1326                IF ( skip_time_do2d_xy /= 0.0_wp )  THEN
[1]1327                   WRITE ( io, 339 )  skip_time_do2d_xy
1328                ENDIF
1329             ELSE
1330                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
1331                                   TRIM( begin_chr ), averaging_interval, &
1332                                   dt_averaging_input, 'k', TRIM( slices ), &
1333                                   TRIM( coordinates )
[1353]1334                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1335                   WRITE ( io, 339 )  skip_time_data_output_av
1336                ENDIF
1337             ENDIF
[1308]1338             IF ( netcdf_data_format > 4 )  THEN
1339                WRITE ( io, 352 )  ntdim_2d_xy(av)
1340             ELSE
1341                WRITE ( io, 353 )
1342             ENDIF
[1]1343          ENDIF
1344
1345          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
1346             i = 1
1347             slices = '/'
1348             coordinates = '/'
1349!
[1551]1350!--          Building strings with index and coordinate information of the
[1]1351!--          slices
1352             DO  WHILE ( section(i,2) /= -9999 )
1353
1354                WRITE (section_chr,'(I5)')  section(i,2)
1355                section_chr = ADJUSTL( section_chr )
1356                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1357
1358                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
1359                coor_chr = ADJUSTL( coor_chr )
1360                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1361
1362                i = i + 1
1363             ENDDO
1364             IF ( av == 0 )  THEN
1365                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
1366                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
1367                                   TRIM( coordinates )
[1353]1368                IF ( skip_time_do2d_xz /= 0.0_wp )  THEN
[1]1369                   WRITE ( io, 339 )  skip_time_do2d_xz
1370                ENDIF
1371             ELSE
1372                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
1373                                   TRIM( begin_chr ), averaging_interval, &
1374                                   dt_averaging_input, 'j', TRIM( slices ), &
1375                                   TRIM( coordinates )
[1353]1376                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1377                   WRITE ( io, 339 )  skip_time_data_output_av
1378                ENDIF
1379             ENDIF
[1308]1380             IF ( netcdf_data_format > 4 )  THEN
1381                WRITE ( io, 352 )  ntdim_2d_xz(av)
1382             ELSE
1383                WRITE ( io, 353 )
1384             ENDIF
[1]1385          ENDIF
1386
1387          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
1388             i = 1
1389             slices = '/'
1390             coordinates = '/'
1391!
[1551]1392!--          Building strings with index and coordinate information of the
[1]1393!--          slices
1394             DO  WHILE ( section(i,3) /= -9999 )
1395
1396                WRITE (section_chr,'(I5)')  section(i,3)
1397                section_chr = ADJUSTL( section_chr )
1398                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1399
1400                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
1401                coor_chr = ADJUSTL( coor_chr )
1402                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1403
1404                i = i + 1
1405             ENDDO
1406             IF ( av == 0 )  THEN
1407                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
1408                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
1409                                   TRIM( coordinates )
[1353]1410                IF ( skip_time_do2d_yz /= 0.0_wp )  THEN
[1]1411                   WRITE ( io, 339 )  skip_time_do2d_yz
1412                ENDIF
1413             ELSE
1414                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
1415                                   TRIM( begin_chr ), averaging_interval, &
1416                                   dt_averaging_input, 'i', TRIM( slices ), &
1417                                   TRIM( coordinates )
[1353]1418                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1419                   WRITE ( io, 339 )  skip_time_data_output_av
1420                ENDIF
1421             ENDIF
[1308]1422             IF ( netcdf_data_format > 4 )  THEN
1423                WRITE ( io, 352 )  ntdim_2d_yz(av)
1424             ELSE
1425                WRITE ( io, 353 )
1426             ENDIF
[1]1427          ENDIF
1428
1429       ENDIF
1430
1431    ENDDO
1432
1433!
1434!-- 3d-arrays
1435    DO  av = 0, 1
1436
1437       i = 1
1438       do3d_chr = ''
1439       DO  WHILE ( do3d(av,i) /= ' ' )
1440
1441          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
1442          i = i + 1
1443
1444       ENDDO
1445
1446       IF ( do3d_chr /= '' )  THEN
1447          IF ( av == 0 )  THEN
1448             WRITE ( io, 336 )  ''
1449          ELSE
1450             WRITE ( io, 336 )  '(time-averaged)'
1451          ENDIF
1452
[1783]1453          output_format = netcdf_data_format_string
1454          IF ( netcdf_deflate == 0 )  THEN
1455             WRITE ( io, 344 )  output_format
1456          ELSE
1457             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1458          ENDIF
[1]1459
1460          IF ( do3d_at_begin )  THEN
1461             begin_chr = 'and at the start'
1462          ELSE
1463             begin_chr = ''
1464          ENDIF
1465          IF ( av == 0 )  THEN
1466             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
1467                                zu(nz_do3d), nz_do3d
1468          ELSE
1469             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
1470                                TRIM( begin_chr ), averaging_interval, &
1471                                dt_averaging_input, zu(nz_do3d), nz_do3d
1472          ENDIF
1473
[1308]1474          IF ( netcdf_data_format > 4 )  THEN
1475             WRITE ( io, 352 )  ntdim_3d(av)
1476          ELSE
1477             WRITE ( io, 353 )
1478          ENDIF
1479
[1]1480          IF ( av == 0 )  THEN
[1353]1481             IF ( skip_time_do3d /= 0.0_wp )  THEN
[1]1482                WRITE ( io, 339 )  skip_time_do3d
1483             ENDIF
1484          ELSE
[1353]1485             IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1486                WRITE ( io, 339 )  skip_time_data_output_av
1487             ENDIF
1488          ENDIF
1489
1490       ENDIF
1491
1492    ENDDO
1493
1494!
[410]1495!-- masked arrays
1496    IF ( masks > 0 )  WRITE ( io, 345 )  &
1497         mask_scale_x, mask_scale_y, mask_scale_z
1498    DO  mid = 1, masks
1499       DO  av = 0, 1
1500
1501          i = 1
1502          domask_chr = ''
1503          DO  WHILE ( domask(mid,av,i) /= ' ' )
1504             domask_chr = TRIM( domask_chr ) // ' ' //  &
1505                          TRIM( domask(mid,av,i) ) // ','
1506             i = i + 1
1507          ENDDO
1508
1509          IF ( domask_chr /= '' )  THEN
1510             IF ( av == 0 )  THEN
1511                WRITE ( io, 346 )  '', mid
1512             ELSE
1513                WRITE ( io, 346 )  ' (time-averaged)', mid
1514             ENDIF
1515
[1783]1516             output_format = netcdf_data_format_string
[1308]1517!--          Parallel output not implemented for mask data, hence
1518!--          output_format must be adjusted.
1519             IF ( netcdf_data_format == 5 ) output_format = 'netCDF4/HDF5'
1520             IF ( netcdf_data_format == 6 ) output_format = 'netCDF4/HDF5 classic'
[1783]1521             IF ( netcdf_deflate == 0 )  THEN
1522                WRITE ( io, 344 )  output_format
1523             ELSE
1524                WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1525             ENDIF
[410]1526
1527             IF ( av == 0 )  THEN
1528                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1529             ELSE
1530                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1531                                   averaging_interval, dt_averaging_input
1532             ENDIF
1533
1534             IF ( av == 0 )  THEN
[1353]1535                IF ( skip_time_domask(mid) /= 0.0_wp )  THEN
[410]1536                   WRITE ( io, 339 )  skip_time_domask(mid)
1537                ENDIF
1538             ELSE
[1353]1539                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[410]1540                   WRITE ( io, 339 )  skip_time_data_output_av
1541                ENDIF
1542             ENDIF
1543!
1544!--          output locations
1545             DO  dim = 1, 3
[1353]1546                IF ( mask(mid,dim,1) >= 0.0_wp )  THEN
[410]1547                   count = 0
[1353]1548                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0_wp )
[410]1549                      count = count + 1
1550                   ENDDO
1551                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1552                                      mask(mid,dim,:count)
[1353]1553                ELSEIF ( mask_loop(mid,dim,1) < 0.0_wp .AND.  &
1554                         mask_loop(mid,dim,2) < 0.0_wp .AND.  &
1555                         mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1556                   WRITE ( io, 350 )  dir(dim), dir(dim)
[1353]1557                ELSEIF ( mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1558                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1559                                      mask_loop(mid,dim,1:2)
1560                ELSE
1561                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1562                                      mask_loop(mid,dim,1:3)
1563                ENDIF
1564             ENDDO
1565          ENDIF
1566
1567       ENDDO
1568    ENDDO
1569
1570!
[1]1571!-- Timeseries
[1322]1572    IF ( dt_dots /= 9999999.9_wp )  THEN
[1]1573       WRITE ( io, 340 )
1574
[1783]1575       output_format = netcdf_data_format_string
1576       IF ( netcdf_deflate == 0 )  THEN
1577          WRITE ( io, 344 )  output_format
1578       ELSE
1579          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1580       ENDIF
[1]1581       WRITE ( io, 341 )  dt_dots
1582    ENDIF
1583
1584#if defined( __dvrp_graphics )
1585!
1586!-- Dvrp-output
[1322]1587    IF ( dt_dvrp /= 9999999.9_wp )  THEN
[1]1588       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1589                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1590       i = 1
1591       l = 0
[336]1592       m = 0
[1]1593       DO WHILE ( mode_dvrp(i) /= ' ' )
1594          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
[130]1595             READ ( mode_dvrp(i), '(10X,I2)' )  j
[1]1596             l = l + 1
1597             IF ( do3d(0,j) /= ' ' )  THEN
[336]1598                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1599                                   isosurface_color(:,l)
[1]1600             ENDIF
1601          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
[130]1602             READ ( mode_dvrp(i), '(6X,I2)' )  j
[336]1603             m = m + 1
1604             IF ( do2d(0,j) /= ' ' )  THEN
1605                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1606                                   slicer_range_limits_dvrp(:,m)
1607             ENDIF
[1]1608          ENDIF
1609          i = i + 1
1610       ENDDO
[237]1611
[336]1612       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1613                          superelevation_y, superelevation, clip_dvrp_l, &
1614                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1615
1616       IF ( TRIM( topography ) /= 'flat' )  THEN
1617          WRITE ( io, 366 )  topography_color
1618          IF ( cluster_size > 1 )  THEN
1619             WRITE ( io, 367 )  cluster_size
1620          ENDIF
[237]1621       ENDIF
1622
[1]1623    ENDIF
1624#endif
1625
1626!
1627!-- Spectra output
[1322]1628    IF ( dt_dosp /= 9999999.9_wp )  THEN
[1]1629       WRITE ( io, 370 )
1630
[1783]1631       output_format = netcdf_data_format_string
1632       IF ( netcdf_deflate == 0 )  THEN
1633          WRITE ( io, 344 )  output_format
1634       ELSE
1635          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1636       ENDIF
[1]1637       WRITE ( io, 371 )  dt_dosp
[1353]1638       IF ( skip_time_dosp /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dosp
[1]1639       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
1640                          ( spectra_direction(i), i = 1,10 ),  &
[189]1641                          ( comp_spectra_level(i), i = 1,100 ), &
1642                          ( plot_spectra_level(i), i = 1,100 ), &
[1]1643                          averaging_interval_sp, dt_averaging_input_pr
1644    ENDIF
1645
1646    WRITE ( io, 99 )
1647
1648!
1649!-- Physical quantities
1650    WRITE ( io, 400 )
1651
1652!
1653!-- Geostrophic parameters
[1551]1654    WRITE ( io, 410 )  phi, omega, f, fs
[1]1655
1656!
1657!-- Other quantities
1658    WRITE ( io, 411 )  g
[1551]1659
[1179]1660    WRITE ( io, 412 )  TRIM( reference_state )
1661    IF ( use_single_reference_value )  THEN
[97]1662       IF ( ocean )  THEN
[1179]1663          WRITE ( io, 413 )  prho_reference
[97]1664       ELSE
[1179]1665          WRITE ( io, 414 )  pt_reference
[97]1666       ENDIF
1667    ENDIF
[1]1668
1669!
1670!-- Cloud physics parameters
[1299]1671    IF ( cloud_physics )  THEN
[57]1672       WRITE ( io, 415 )
1673       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
[1822]1674       IF ( microphysics_seifert )  THEN
[1353]1675          WRITE ( io, 510 ) 1.0E-6_wp * nc_const
[1822]1676          WRITE ( io, 511 ) c_sedimentation
[1115]1677       ENDIF
[1]1678    ENDIF
1679
1680!
[824]1681!-- Cloud physcis parameters / quantities / numerical methods
1682    WRITE ( io, 430 )
1683    IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets)  THEN
1684       WRITE ( io, 431 )
1685    ELSEIF ( humidity  .AND.  cloud_physics )  THEN
1686       WRITE ( io, 432 )
[1496]1687       IF ( cloud_top_radiation )  WRITE ( io, 132 )
[1822]1688       IF ( microphysics_kessler )  THEN
1689          WRITE ( io, 133 )
1690       ELSEIF ( microphysics_seifert )  THEN
[1115]1691          IF ( drizzle )  WRITE ( io, 506 )
[1822]1692          WRITE ( io, 505 )
1693          IF ( turbulence )  WRITE ( io, 507 )
1694          IF ( ventilation_effect )  WRITE ( io, 508 )
1695          IF ( limiter_sedimentation )  WRITE ( io, 509 )
[1115]1696       ENDIF
[824]1697    ELSEIF ( humidity  .AND.  cloud_droplets )  THEN
1698       WRITE ( io, 433 )
1699       IF ( curvature_solution_effects )  WRITE ( io, 434 )
[825]1700       IF ( collision_kernel /= 'none' )  THEN
1701          WRITE ( io, 435 )  TRIM( collision_kernel )
[828]1702          IF ( collision_kernel(6:9) == 'fast' )  THEN
1703             WRITE ( io, 436 )  radius_classes, dissipation_classes
1704          ENDIF
[825]1705       ELSE
[828]1706          WRITE ( io, 437 )
[825]1707       ENDIF
[824]1708    ENDIF
1709
1710!
[1]1711!-- LES / turbulence parameters
1712    WRITE ( io, 450 )
1713
1714!--
1715! ... LES-constants used must still be added here
1716!--
1717    IF ( constant_diffusion )  THEN
1718       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1719                          prandtl_number
1720    ENDIF
1721    IF ( .NOT. constant_diffusion)  THEN
[1353]1722       IF ( e_init > 0.0_wp )  WRITE ( io, 455 )  e_init
1723       IF ( e_min > 0.0_wp )  WRITE ( io, 454 )  e_min
[1]1724       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1725    ENDIF
1726
1727!
1728!-- Special actions during the run
1729    WRITE ( io, 470 )
1730    IF ( create_disturbances )  THEN
1731       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1732                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1733                          zu(disturbance_level_ind_t), disturbance_level_ind_t
[707]1734       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1]1735          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1736       ELSE
1737          WRITE ( io, 473 )  disturbance_energy_limit
1738       ENDIF
1739       WRITE ( io, 474 )  TRIM( random_generator )
1740    ENDIF
[1353]1741    IF ( pt_surface_initial_change /= 0.0_wp )  THEN
[1]1742       WRITE ( io, 475 )  pt_surface_initial_change
1743    ENDIF
[1353]1744    IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1745       WRITE ( io, 476 )  q_surface_initial_change       
1746    ENDIF
[1353]1747    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1748       WRITE ( io, 477 )  q_surface_initial_change       
1749    ENDIF
1750
[60]1751    IF ( particle_advection )  THEN
[1]1752!
[60]1753!--    Particle attributes
1754       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1755                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
[1359]1756                          end_time_prel
[60]1757       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1758       IF ( random_start_position )  WRITE ( io, 481 )
[1575]1759       IF ( seed_follows_topography )  WRITE ( io, 496 )
[60]1760       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1761       WRITE ( io, 495 )  total_number_of_particles
[1322]1762       IF ( dt_write_particle_data /= 9999999.9_wp )  THEN
[60]1763          WRITE ( io, 485 )  dt_write_particle_data
[1327]1764          IF ( netcdf_data_format > 1 )  THEN
1765             output_format = 'netcdf (64 bit offset) and binary'
[1]1766          ELSE
[1327]1767             output_format = 'netcdf and binary'
[1]1768          ENDIF
[1783]1769          IF ( netcdf_deflate == 0 )  THEN
1770             WRITE ( io, 344 )  output_format
1771          ELSE
1772             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1773          ENDIF
[1]1774       ENDIF
[1322]1775       IF ( dt_dopts /= 9999999.9_wp )  WRITE ( io, 494 )  dt_dopts
[60]1776       IF ( write_particle_statistics )  WRITE ( io, 486 )
[1]1777
[60]1778       WRITE ( io, 487 )  number_of_particle_groups
[1]1779
[60]1780       DO  i = 1, number_of_particle_groups
[1322]1781          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9_wp )  THEN
[1353]1782             WRITE ( io, 490 )  i, 0.0_wp
[60]1783             WRITE ( io, 492 )
[1]1784          ELSE
[60]1785             WRITE ( io, 490 )  i, radius(i)
[1353]1786             IF ( density_ratio(i) /= 0.0_wp )  THEN
[60]1787                WRITE ( io, 491 )  density_ratio(i)
1788             ELSE
1789                WRITE ( io, 492 )
1790             ENDIF
[1]1791          ENDIF
[60]1792          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1793                             pdx(i), pdy(i), pdz(i)
[336]1794          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
[60]1795       ENDDO
[1]1796
[60]1797    ENDIF
[1]1798
[60]1799
[1]1800!
1801!-- Parameters of 1D-model
1802    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1803       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1804                          mixing_length_1d, dissipation_1d
1805       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1806          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1807       ENDIF
1808    ENDIF
1809
1810!
[1551]1811!-- User-defined information
[1]1812    CALL user_header( io )
1813
1814    WRITE ( io, 99 )
1815
1816!
1817!-- Write buffer contents to disc immediately
[1808]1818    FLUSH( io )
[1]1819
1820!
1821!-- Here the FORMATs start
1822
1823 99 FORMAT (1X,78('-'))
[1468]1824100 FORMAT (/1X,'******************************',4X,44('-')/        &
1825            1X,'* ',A,' *',4X,A/                               &
1826            1X,'******************************',4X,44('-'))
1827101 FORMAT (35X,'coupled run using MPI-',I1,': ',A/ &
1828            35X,42('-'))
1829102 FORMAT (/' Date:                 ',A8,4X,'Run:       ',A20/      &
1830            ' Time:                 ',A8,4X,'Run-No.:   ',I2.2/     &
[1106]1831            ' Run on host:        ',A10)
[1]1832#if defined( __parallel )
[1468]1833103 FORMAT (' Number of PEs:',10X,I6,4X,'Processor grid (x,y): (',I4,',',I4, &
[1]1834              ')',1X,A)
[1468]1835104 FORMAT (' Number of PEs:',10X,I6,4X,'Tasks:',I4,'   threads per task:',I4/ &
1836              35X,'Processor grid (x,y): (',I4,',',I4,')',1X,A)
1837105 FORMAT (35X,'One additional PE is used to handle'/37X,'the dvrp output!')
1838106 FORMAT (35X,'A 1d-decomposition along x is forced'/ &
1839            35X,'because the job is running on an SMP-cluster')
1840107 FORMAT (35X,'A 1d-decomposition along ',A,' is used')
1841108 FORMAT (35X,'Max. # of parallel I/O streams is ',I5)
1842109 FORMAT (35X,'Precursor run for coupled atmos-ocean run'/ &
1843            35X,42('-'))
1844114 FORMAT (35X,'Coupled atmosphere-ocean run following'/ &
1845            35X,'independent precursor runs'/             &
1846            35X,42('-'))
[1111]1847117 FORMAT (' Accelerator boards / node:  ',I2)
[1]1848#endif
1849110 FORMAT (/' Numerical Schemes:'/ &
1850             ' -----------------'/)
1851111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1852112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
[1697]1853            '     Iterations (initial/other): ',I3,'/',I3,'  omega =',F6.3)
[1]1854113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1855                  ' or Upstream')
[1216]1856115 FORMAT ('     FFT and transpositions are overlapping')
[1]1857116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1858                  ' or Upstream')
1859118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
[1106]1860119 FORMAT (' --> Galilei-Transform applied to horizontal advection:'/ &
1861            '     translation velocity = ',A/ &
[1]1862            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
[1111]1863120 FORMAT (' Accelerator boards: ',8X,I2)
[1]1864122 FORMAT (' --> Time differencing scheme: ',A)
[108]1865123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
[1697]1866            '     maximum damping coefficient:',F6.3, ' 1/s')
[1]1867129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1868130 FORMAT (' --> Additional prognostic equation for the total water content')
[940]1869131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
1870                  F6.2, ' K assumed')
[824]1871132 FORMAT ('     Parameterization of long-wave radiation processes via'/ &
[1]1872            '     effective emissivity scheme')
[824]1873133 FORMAT ('     Precipitation parameterization via Kessler-Scheme')
[1]1874134 FORMAT (' --> Additional prognostic equation for a passive scalar')
[1575]1875135 FORMAT (' --> Solve perturbation pressure via ',A,' method (', &
[1]1876                  A,'-cycle)'/ &
1877            '     number of grid levels:                   ',I2/ &
1878            '     Gauss-Seidel red/black iterations:       ',I2)
1879136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1880                  I3,')')
1881137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1882            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1883                  I3,')'/ &
1884            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1885                  I3,')')
[63]1886139 FORMAT (' --> Loop optimization method: ',A)
[1]1887140 FORMAT ('     maximum residual allowed:                ',E10.3)
1888141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1889142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1890                  'step')
[87]1891143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1892                  'kinetic energy')
[927]1893144 FORMAT ('     masking method is used')
[1]1894150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
[241]1895                  'conserved'/ &
1896            '     using the ',A,' mode')
1897151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
[306]1898152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
1899           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
1900           /'     starting from dp_level_b =', F8.3, 'm', A /)
[1365]1901160 FORMAT (//' Large scale forcing and nudging:'/ &
1902              ' -------------------------------'/)
1903161 FORMAT (' --> No large scale forcing from external is used (default) ')
1904162 FORMAT (' --> Large scale forcing from external file LSF_DATA is used: ')
1905163 FORMAT ('     - large scale advection tendencies ')
1906164 FORMAT ('     - large scale subsidence velocity w_subs ')
1907165 FORMAT ('     - large scale subsidence tendencies ')
1908167 FORMAT ('     - and geostrophic wind components ug and vg')
1909168 FORMAT (' --> Large-scale vertical motion is used in the ', &
[1299]1910                  'prognostic equation(s) for')
[1365]1911169 FORMAT ('     the scalar(s) only')
1912170 FORMAT (' --> Nudging is used')
1913171 FORMAT (' --> No nudging is used (default) ')
1914180 FORMAT ('     - prescribed surface values for temperature')
[1376]1915181 FORMAT ('     - prescribed surface fluxes for temperature')
1916182 FORMAT ('     - prescribed surface values for humidity')
[1365]1917183 FORMAT ('     - prescribed surface fluxes for humidity')
[1]1918200 FORMAT (//' Run time and time step information:'/ &
1919             ' ----------------------------------'/)
[1106]1920201 FORMAT ( ' Timestep:             variable     maximum value: ',F6.3,' s', &
[1697]1921             '    CFL-factor:',F5.2)
[1106]1922202 FORMAT ( ' Timestep:          dt = ',F6.3,' s'/)
1923203 FORMAT ( ' Start time:          ',F9.3,' s'/ &
1924             ' End time:            ',F9.3,' s')
[1]1925204 FORMAT ( A,F9.3,' s')
1926205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
[1106]1927206 FORMAT (/' Time reached:        ',F9.3,' s'/ &
1928             ' CPU-time used:       ',F9.3,' s     per timestep:               ', &
1929               '  ',F9.3,' s'/                                                    &
[1111]1930             '                                      per second of simulated tim', &
[1]1931               'e: ',F9.3,' s')
[1106]1932207 FORMAT ( ' Coupling start time: ',F9.3,' s')
[1]1933250 FORMAT (//' Computational grid and domain size:'/ &
1934              ' ----------------------------------'// &
1935              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
1936              ' m    dz =    ',F7.3,' m'/ &
1937              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1938              ' m  z(u) = ',F10.3,' m'/)
1939252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
[1697]1940              ' factor:',F6.3/ &
[1]1941            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
1942254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1943            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1944260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1945             ' degrees')
[1551]1946270 FORMAT (//' Topography information:'/ &
1947              ' ----------------------'// &
[1]1948              1X,'Topography: ',A)
1949271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1950              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1951                ' / ',I4)
[240]1952272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
1953              ' direction' / &
1954              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
1955              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
[256]1956278 FORMAT (' Topography grid definition convention:'/ &
1957            ' cell edge (staggered grid points'/  &
1958            ' (u in x-direction, v in y-direction))' /)
1959279 FORMAT (' Topography grid definition convention:'/ &
1960            ' cell center (scalar grid points)' /)
[1]1961300 FORMAT (//' Boundary conditions:'/ &
1962             ' -------------------'// &
1963             '                     p                    uv             ', &
[1551]1964             '                     pt'// &
[1]1965             ' B. bound.: ',A/ &
1966             ' T. bound.: ',A)
[97]1967301 FORMAT (/'                     ',A// &
[1]1968             ' B. bound.: ',A/ &
1969             ' T. bound.: ',A)
[19]1970303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1971304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1972305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
1973               'computational u,v-level:'// &
[1697]1974             '       zp = ',F6.2,' m   z0 =',F7.4,' m   z0h =',F8.5,&
1975             ' m   kappa =',F5.2/ &
1976             '       Rif value range:   ',F8.2,' <= rif <=',F6.2)
[97]1977306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
[1]1978307 FORMAT ('       Heatflux has a random normal distribution')
1979308 FORMAT ('       Predefined surface temperature')
[97]1980309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
[1]1981310 FORMAT (//'    1D-Model:'// &
1982             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1983311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
1984312 FORMAT ('       Predefined surface humidity')
1985313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1986314 FORMAT ('       Predefined scalar value at the surface')
[19]1987315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
[102]1988316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1989                    'atmosphere model')
[1]1990317 FORMAT (//' Lateral boundaries:'/ &
1991            '       left/right:  ',A/    &
1992            '       north/south: ',A)
[1159]1993318 FORMAT (/'       use_cmax: ',L1 / &
1994            '       pt damping layer width = ',F8.2,' m, pt ', &
[1697]1995                    'damping factor =',F7.4)
[151]1996319 FORMAT ('       turbulence recycling at inflow switched on'/ &
1997            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1998            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
1999320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
[103]2000            '                                          v: ',F9.6,' m**2/s**2')
[1365]2001321 FORMAT (//' Initial profiles:'/ &
2002              ' ----------------')
[1560]2003322 FORMAT ('       turbulence recycling at inflow switched on'/ &
2004            '       y shift of the recycled inflow turbulence switched on'/ &
2005            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
[1592]2006            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m'/)
[151]2007325 FORMAT (//' List output:'/ &
[1]2008             ' -----------'//  &
2009            '    1D-Profiles:'/    &
2010            '       Output every             ',F8.2,' s')
[151]2011326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
[1]2012            '       Averaging input every    ',F8.2,' s')
2013330 FORMAT (//' Data output:'/ &
2014             ' -----------'/)
2015331 FORMAT (/'    1D-Profiles:')
2016332 FORMAT (/'       ',A)
2017333 FORMAT ('       Output every             ',F8.2,' s',/ &
2018            '       Time averaged over       ',F8.2,' s'/ &
2019            '       Averaging input every    ',F8.2,' s')
2020334 FORMAT (/'    2D-Arrays',A,':')
2021335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2022            '       Output every             ',F8.2,' s  ',A/ &
2023            '       Cross sections at ',A1,' = ',A/ &
2024            '       scalar-coordinates:   ',A,' m'/)
2025336 FORMAT (/'    3D-Arrays',A,':')
2026337 FORMAT (/'       Arrays: ',A/ &
2027            '       Output every             ',F8.2,' s  ',A/ &
2028            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
2029339 FORMAT ('       No output during initial ',F8.2,' s')
2030340 FORMAT (/'    Time series:')
2031341 FORMAT ('       Output every             ',F8.2,' s'/)
2032342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2033            '       Output every             ',F8.2,' s  ',A/ &
2034            '       Time averaged over       ',F8.2,' s'/ &
2035            '       Averaging input every    ',F8.2,' s'/ &
2036            '       Cross sections at ',A1,' = ',A/ &
2037            '       scalar-coordinates:   ',A,' m'/)
2038343 FORMAT (/'       Arrays: ',A/ &
2039            '       Output every             ',F8.2,' s  ',A/ &
2040            '       Time averaged over       ',F8.2,' s'/ &
2041            '       Averaging input every    ',F8.2,' s'/ &
2042            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
[292]2043344 FORMAT ('       Output format: ',A/)
[410]2044345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
2045            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
2046            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
2047            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
2048346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
2049347 FORMAT ('       Variables: ',A/ &
2050            '       Output every             ',F8.2,' s')
2051348 FORMAT ('       Variables: ',A/ &
2052            '       Output every             ',F8.2,' s'/ &
2053            '       Time averaged over       ',F8.2,' s'/ &
2054            '       Averaging input every    ',F8.2,' s')
2055349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2056            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
2057            13('       ',8(F8.2,',')/) )
2058350 FORMAT (/'       Output locations in ',A,'-direction: ', &
2059            'all gridpoints along ',A,'-direction (default).' )
2060351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2061            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
2062            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
[1313]2063352 FORMAT  (/'       Number of output time levels allowed: ',I3 /)
2064353 FORMAT  (/'       Number of output time levels allowed: unlimited' /)
[1783]2065354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
[1]2066#if defined( __dvrp_graphics )
2067360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
2068            '       Output every      ',F7.1,' s'/ &
2069            '       Output mode:      ',A/ &
2070            '       Host / User:      ',A,' / ',A/ &
2071            '       Directory:        ',A// &
2072            '       The sequence contains:')
[337]2073361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
2074            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
2075362 FORMAT (/'       Slicer plane ',A/ &
[336]2076            '       Slicer limits: [',F6.2,',',F6.2,']')
[337]2077365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
[336]2078            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
2079                     ')'/ &
2080            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
2081            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
[337]2082366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
[336]2083367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
[1]2084#endif
2085370 FORMAT ('    Spectra:')
2086371 FORMAT ('       Output every ',F7.1,' s'/)
2087372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
2088            '       Directions: ', 10(A5,',')/                         &
[189]2089            '       height levels  k = ', 20(I3,',')/                  &
2090            '                          ', 20(I3,',')/                  &
2091            '                          ', 20(I3,',')/                  &
2092            '                          ', 20(I3,',')/                  &
2093            '                          ', 19(I3,','),I3,'.'/           &
[1]2094            '       height levels selected for standard plot:'/        &
[189]2095            '                      k = ', 20(I3,',')/                  &
2096            '                          ', 20(I3,',')/                  &
2097            '                          ', 20(I3,',')/                  &
2098            '                          ', 20(I3,',')/                  &
2099            '                          ', 19(I3,','),I3,'.'/           &
[1]2100            '       Time averaged over ', F7.1, ' s,' /                &
2101            '       Profiles for the time averaging are taken every ', &
2102                    F6.1,' s')
2103400 FORMAT (//' Physical quantities:'/ &
2104              ' -------------------'/)
[1551]2105410 FORMAT ('    Geograph. latitude  :   phi    = ',F4.1,' degr'/   &
[1697]2106            '    Angular velocity    :   omega  =',E10.3,' rad/s'/  &
[1551]2107            '    Coriolis parameter  :   f      = ',F9.6,' 1/s'/    &
2108            '                            f*     = ',F9.6,' 1/s')
2109411 FORMAT (/'    Gravity             :   g      = ',F4.1,' m/s**2')
[1179]2110412 FORMAT (/'    Reference state used in buoyancy terms: ',A)
2111413 FORMAT ('       Reference density in buoyancy terms: ',F8.3,' kg/m**3')
2112414 FORMAT ('       Reference temperature in buoyancy terms: ',F8.4,' K')
[1551]2113415 FORMAT (/' Cloud physics parameters:'/ &
2114             ' ------------------------'/)
2115416 FORMAT ('    Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
2116            '    Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
[1697]2117            '    Density of air     :   rho_0 =',F6.3,' kg/m**3'/  &
[1551]2118            '    Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
[1697]2119            '    Vapourization heat :   L_v   =',E9.2,' J/kg')
[1551]2120417 FORMAT ('    Geograph. longitude :   lambda = ',F4.1,' degr')
2121418 FORMAT (/'    Day of the year at model start :   day_init      =     ',I3 &
2122            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
[1]2123420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
2124            '       Height:        ',A,'  m'/ &
2125            '       Temperature:   ',A,'  K'/ &
2126            '       Gradient:      ',A,'  K/100m'/ &
2127            '       Gridpoint:     ',A)
2128421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
2129            '       Height:      ',A,'  m'/ &
2130            '       Humidity:    ',A,'  kg/kg'/ &
2131            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
2132            '       Gridpoint:   ',A)
2133422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
2134            '       Height:                  ',A,'  m'/ &
2135            '       Scalar concentration:    ',A,'  kg/m**3'/ &
2136            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
2137            '       Gridpoint:               ',A)
2138423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
2139            '       Height:      ',A,'  m'/ &
2140            '       ug:          ',A,'  m/s'/ &
2141            '       Gradient:    ',A,'  1/100s'/ &
2142            '       Gridpoint:   ',A)
2143424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
2144            '       Height:      ',A,'  m'/ &
[97]2145            '       vg:          ',A,'  m/s'/ &
[1]2146            '       Gradient:    ',A,'  1/100s'/ &
2147            '       Gridpoint:   ',A)
[97]2148425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
2149            '       Height:     ',A,'  m'/ &
2150            '       Salinity:   ',A,'  psu'/ &
2151            '       Gradient:   ',A,'  psu/100m'/ &
2152            '       Gridpoint:  ',A)
[411]2153426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
2154            '       Height:      ',A,'  m'/ &
2155            '       w_subs:      ',A,'  m/s'/ &
2156            '       Gradient:    ',A,'  (m/s)/100m'/ &
2157            '       Gridpoint:   ',A)
[767]2158427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
2159                  ' profiles')
[1241]2160428 FORMAT (/'    Initial profiles (u, v, pt, q) are taken from file '/ &
2161             '    NUDGING_DATA')
[824]2162430 FORMAT (//' Cloud physics quantities / methods:'/ &
2163              ' ----------------------------------'/)
2164431 FORMAT ('    Humidity is treated as purely passive scalar (no condensati', &
2165                 'on)')
2166432 FORMAT ('    Bulk scheme with liquid water potential temperature and'/ &
2167            '    total water content is used.'/ &
2168            '    Condensation is parameterized via 0% - or 100% scheme.')
2169433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
2170                 'icle model')
2171434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
2172                 ' droplets < 1.0E-6 m')
[825]2173435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
[828]2174436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
2175                    'are used'/ &
2176            '          number of radius classes:       ',I3,'    interval ', &
2177                       '[1.0E-6,2.0E-4] m'/ &
2178            '          number of dissipation classes:   ',I2,'    interval ', &
2179                       '[0,1000] cm**2/s**3')
2180437 FORMAT ('    Droplet collision is switched off')
[1]2181450 FORMAT (//' LES / Turbulence quantities:'/ &
2182              ' ---------------------------'/)
[824]2183451 FORMAT ('    Diffusion coefficients are constant:'/ &
2184            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
[1697]2185453 FORMAT ('    Mixing length is limited to',F5.2,' * z')
[824]2186454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
2187455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
[1]2188470 FORMAT (//' Actions during the simulation:'/ &
2189              ' -----------------------------'/)
[94]2190471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
[1697]2191            '    Disturbance amplitude           :    ',F5.2, ' m/s'/       &
[94]2192            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
2193            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
[1]2194472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
2195                 ' to i/j =',I4)
2196473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
[1697]2197                 F6.3, ' m**2/s**2')
[1]2198474 FORMAT ('    Random number generator used    : ',A/)
2199475 FORMAT ('    The surface temperature is increased (or decreased, ', &
2200                 'respectively, if'/ &
2201            '    the value is negative) by ',F5.2,' K at the beginning of the',&
2202                 ' 3D-simulation'/)
2203476 FORMAT ('    The surface humidity is increased (or decreased, ',&
2204                 'respectively, if the'/ &
2205            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
2206                 ' the 3D-simulation'/)
2207477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
2208                 'respectively, if the'/ &
2209            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
2210                 ' the 3D-simulation'/)
2211480 FORMAT ('    Particles:'/ &
2212            '    ---------'// &
2213            '       Particle advection is active (switched on at t = ', F7.1, &
2214                    ' s)'/ &
2215            '       Start of new particle generations every  ',F6.1,' s'/ &
2216            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
2217            '                            bottom:     ', A, ' top:         ', A/&
2218            '       Maximum particle age:                 ',F9.1,' s'/ &
[1359]2219            '       Advection stopped at t = ',F9.1,' s'/)
[1]2220481 FORMAT ('       Particles have random start positions'/)
[336]2221482 FORMAT ('          Particles are advected only horizontally'/)
[1]2222485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
2223486 FORMAT ('       Particle statistics are written on file'/)
2224487 FORMAT ('       Number of particle groups: ',I2/)
2225488 FORMAT ('       SGS velocity components are used for particle advection'/ &
[1697]2226            '          minimum timestep for advection:', F8.5/)
[1]2227489 FORMAT ('       Number of particles simultaneously released at each ', &
2228                    'point: ', I5/)
2229490 FORMAT ('       Particle group ',I2,':'/ &
2230            '          Particle radius: ',E10.3, 'm')
2231491 FORMAT ('          Particle inertia is activated'/ &
[1697]2232            '             density_ratio (rho_fluid/rho_particle) =',F6.3/)
[1]2233492 FORMAT ('          Particles are advected only passively (no inertia)'/)
2234493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
2235            '                                         y:',F8.1,' - ',F8.1,' m'/&
2236            '                                         z:',F8.1,' - ',F8.1,' m'/&
2237            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
2238                       ' m  dz = ',F8.1,' m'/)
2239494 FORMAT ('       Output of particle time series in NetCDF format every ', &
2240                    F8.2,' s'/)
2241495 FORMAT ('       Number of particles in total domain: ',I10/)
[1575]2242496 FORMAT ('       Initial vertical particle positions are interpreted ', &
2243                    'as relative to the given topography')
[1]2244500 FORMAT (//' 1D-Model parameters:'/                           &
2245              ' -------------------'//                           &
2246            '    Simulation time:                   ',F8.1,' s'/ &
2247            '    Run-controll output every:         ',F8.1,' s'/ &
2248            '    Vertical profile output every:     ',F8.1,' s'/ &
2249            '    Mixing length calculation:         ',A/         &
2250            '    Dissipation calculation:           ',A/)
2251502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
[667]2252503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2253504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
[1115]2254505 FORMAT ('    Precipitation parameterization via Seifert-Beheng-Scheme')
2255506 FORMAT ('    Drizzle parameterization via Stokes law')
2256507 FORMAT ('    Turbulence effects on precipitation process')
2257508 FORMAT ('    Ventilation effects on evaporation of rain drops')
2258509 FORMAT ('    Slope limiter used for sedimentation process')
[1551]2259510 FORMAT ('    Droplet density    :   N_c   = ',F6.1,' 1/cm**3')
2260511 FORMAT ('    Sedimentation Courant number:                  '/&
[1697]2261            '                               C_s   =',F4.1,'        ')
[1429]2262512 FORMAT (/' Date:                 ',A8,6X,'Run:       ',A20/      &
2263            ' Time:                 ',A8,6X,'Run-No.:   ',I2.2/     &
2264            ' Run on host:        ',A10,6X,'En-No.:    ',I2.2)
[1557]2265513 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order ' // & 
2266            '+ monotonic adjustment')
[1791]2267600 FORMAT (/' Nesting informations:'/ &
2268            ' --------------------'/ &
[1797]2269            ' Nesting mode:                     ',A/ &
2270            ' Nesting-datatransfer mode:        ',A// &
[1791]2271            ' Nest id  parent  number   lower left coordinates   name'/ &
2272            ' (*=me)     id    of PEs      x (m)     y (m)' )
2273601 FORMAT (2X,A1,1X,I2.2,6X,I2.2,5X,I5,5X,F8.2,2X,F8.2,5X,A)
[1]2274
2275 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.