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

Last change on this file since 1791 was 1791, checked in by raasch, 8 years ago

output of nesting informations of all domains; filling up redundant ghost points; renaming of variables, etc.; formatting cleanup

  • Property svn:keywords set to Id
File size: 95.9 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!
[1691]16! Copyright 1997-2015 Leibniz Universitaet Hannover
[1036]17!--------------------------------------------------------------------------------!
18!
[254]19! Current revisions:
[1]20! -----------------
[1791]21! output of nesting informations of all domains
22!
[1485]23! Former revisions:
24! -----------------
25! $Id: header.f90 1791 2016-03-11 10:41:25Z raasch $
26!
[1789]27! 1788 2016-03-10 11:01:04Z maronga
28! Parameter dewfall removed
29!
[1787]30! 1786 2016-03-08 05:49:27Z raasch
31! cpp-direktives for spectra removed
32!
[1784]33! 1783 2016-03-06 18:36:17Z raasch
34! netcdf module and variable names changed, output of netcdf_deflate
35!
[1765]36! 1764 2016-02-28 12:45:19Z raasch
37! output of nesting informations
38!
[1698]39! 1697 2015-10-28 17:14:10Z raasch
40! small E- and F-FORMAT changes to avoid informative compiler messages about
41! insufficient field width
42!
[1692]43! 1691 2015-10-26 16:17:44Z maronga
44! Renamed prandtl_layer to constant_flux_layer, renames rif_min/rif_max to
45! zeta_min/zeta_max.
46!
[1683]47! 1682 2015-10-07 23:56:08Z knoop
48! Code annotations made doxygen readable
49!
[1676]50! 1675 2015-10-02 08:28:59Z gronemeier
51! Bugfix: Definition of topography grid levels
52!
[1662]53! 1660 2015-09-21 08:15:16Z gronemeier
54! Bugfix: Definition of building/street canyon height if vertical grid stretching
55!         starts below the maximum topography height.
56!
[1591]57! 1590 2015-05-08 13:56:27Z maronga
58! Bugfix: Added TRIM statements for character strings for LSM and radiation code
59!
[1586]60! 1585 2015-04-30 07:05:52Z maronga
61! Further output for radiation model(s).
62!
[1576]63! 1575 2015-03-27 09:56:27Z raasch
64! adjustments for psolver-queries, output of seed_follows_topography
65!
[1561]66! 1560 2015-03-06 10:48:54Z keck
67! output for recycling y shift
68!
[1558]69! 1557 2015-03-05 16:43:04Z suehring
70! output for monotonic limiter
71!
[1552]72! 1551 2015-03-03 14:18:16Z maronga
73! Added informal output for land surface model and radiation model. Removed typo.
74!
[1497]75! 1496 2014-12-02 17:25:50Z maronga
76! Renamed: "radiation -> "cloud_top_radiation"
77!
[1485]78! 1484 2014-10-21 10:53:05Z kanani
[1484]79! Changes due to new module structure of the plant canopy model:
80!   module plant_canopy_model_mod and output for new canopy model parameters
81!   (alpha_lad, beta_lad, lai_beta,...) added,
82!   drag_coefficient, leaf_surface_concentration and scalar_exchange_coefficient
83!   renamed to canopy_drag_coeff, leaf_surface_conc and leaf_scalar_exch_coeff,
84!   learde renamed leaf_area_density.
85! Bugfix: DO-WHILE-loop for lad header information additionally restricted
86! by maximum number of gradient levels (currently 10)
[1483]87!
88! 1482 2014-10-18 12:34:45Z raasch
89! information about calculated or predefined virtual processor topology adjusted
90!
[1469]91! 1468 2014-09-24 14:06:57Z maronga
92! Adapted for use on up to 6-digit processor cores
93!
[1430]94! 1429 2014-07-15 12:53:45Z knoop
95! header exended to provide ensemble_member_nr if specified
96!
[1377]97! 1376 2014-04-26 11:21:22Z boeske
98! Correction of typos
99!
[1366]100! 1365 2014-04-22 15:03:56Z boeske
101! New section 'Large scale forcing and nudging':
102! output of large scale forcing and nudging information,
103! new section for initial profiles created
104!
[1360]105! 1359 2014-04-11 17:15:14Z hoffmann
106! dt_sort_particles removed
107!
[1354]108! 1353 2014-04-08 15:21:23Z heinze
109! REAL constants provided with KIND-attribute
110!
[1329]111! 1327 2014-03-21 11:00:16Z raasch
112! parts concerning iso2d and avs output removed,
113! -netcdf output queries
114!
[1325]115! 1324 2014-03-21 09:13:16Z suehring
116! Bugfix: module spectrum added
117!
[1323]118! 1322 2014-03-20 16:38:49Z raasch
119! REAL functions provided with KIND-attribute,
120! some REAL constants defined as wp-kind
121!
[1321]122! 1320 2014-03-20 08:40:49Z raasch
[1320]123! ONLY-attribute added to USE-statements,
124! kind-parameters added to all INTEGER and REAL declaration statements,
125! kinds are defined in new module kinds,
126! revision history before 2012 removed,
127! comment fields (!:) to be used for variable explanations added to
128! all variable declaration statements
[1321]129!
[1309]130! 1308 2014-03-13 14:58:42Z fricke
131! output of the fixed number of output time levels
132! output_format adjusted for masked data if netcdf_data_format > 5
133!
[1300]134! 1299 2014-03-06 13:15:21Z heinze
135! output for using large_scale subsidence in combination
136! with large_scale_forcing
137! reformatting, more detailed explanations
138!
[1242]139! 1241 2013-10-30 11:36:58Z heinze
140! output for nudging + large scale forcing from external file
141!
[1217]142! 1216 2013-08-26 09:31:42Z raasch
143! output for transpose_compute_overlap
144!
[1213]145! 1212 2013-08-15 08:46:27Z raasch
146! output for poisfft_hybrid removed
147!
[1182]148! 1179 2013-06-14 05:57:58Z raasch
149! output of reference_state, use_reference renamed use_single_reference_value
150!
[1160]151! 1159 2013-05-21 11:58:22Z fricke
152! +use_cmax
153!
[1116]154! 1115 2013-03-26 18:16:16Z hoffmann
155! descriptions for Seifert-Beheng-cloud-physics-scheme added
156!
[1112]157! 1111 2013-03-08 23:54:10Z raasch
158! output of accelerator board information
159! ibc_p_b = 2 removed
160!
[1109]161! 1108 2013-03-05 07:03:32Z raasch
162! bugfix for r1106
163!
[1107]164! 1106 2013-03-04 05:31:38Z raasch
165! some format changes for coupled runs
166!
[1093]167! 1092 2013-02-02 11:24:22Z raasch
168! unused variables removed
169!
[1037]170! 1036 2012-10-22 13:43:42Z raasch
171! code put under GPL (PALM 3.9)
172!
[1035]173! 1031 2012-10-19 14:35:30Z raasch
174! output of netCDF data format modified
175!
[1017]176! 1015 2012-09-27 09:23:24Z raasch
[1365]177! output of Adjustment of mixing length to the Prandtl mixing length at first
[1017]178! grid point above ground removed
179!
[1004]180! 1003 2012-09-14 14:35:53Z raasch
181! output of information about equal/unequal subdomain size removed
182!
[1002]183! 1001 2012-09-13 14:08:46Z raasch
184! all actions concerning leapfrog- and upstream-spline-scheme removed
185!
[979]186! 978 2012-08-09 08:28:32Z fricke
187! -km_damp_max, outflow_damping_width
188! +pt_damping_factor, pt_damping_width
189! +z0h
190!
[965]191! 964 2012-07-26 09:14:24Z raasch
192! output of profil-related quantities removed
193!
[941]194! 940 2012-07-09 14:31:00Z raasch
195! Output in case of simulations for pure neutral stratification (no pt-equation
196! solved)
197!
[928]198! 927 2012-06-06 19:15:04Z raasch
199! output of masking_method for mg-solver
200!
[869]201! 868 2012-03-28 12:21:07Z raasch
202! translation velocity in Galilean transformation changed to 0.6 * ug
203!
[834]204! 833 2012-02-22 08:55:55Z maronga
205! Adjusted format for leaf area density
206!
[829]207! 828 2012-02-21 12:00:36Z raasch
208! output of dissipation_classes + radius_classes
209!
[826]210! 825 2012-02-19 03:03:44Z raasch
211! Output of cloud physics parameters/quantities complemented and restructured
212!
[1]213! Revision 1.1  1997/08/11 06:17:20  raasch
214! Initial revision
215!
216!
217! Description:
218! ------------
[1764]219!> Writing a header with all important information about the current run.
[1682]220!> This subroutine is called three times, two times at the beginning
221!> (writing information on files RUN_CONTROL and HEADER) and one time at the
222!> end of the run, then writing additional information about CPU-usage on file
223!> header.
[411]224!-----------------------------------------------------------------------------!
[1682]225 SUBROUTINE header
226 
[1]227
[1320]228    USE arrays_3d,                                                             &
[1660]229        ONLY:  pt_init, qsws, q_init, sa_init, shf, ug, vg, w_subs, zu, zw
[1320]230       
[1]231    USE control_parameters
[1320]232       
233    USE cloud_parameters,                                                      &
234        ONLY:  cp, curvature_solution_effects, c_sedimentation,                &
235               limiter_sedimentation, l_v, nc_const, r_d, ventilation_effect
236       
237    USE cpulog,                                                                &
238        ONLY:  log_point_s
239       
240    USE dvrp_variables,                                                        &
241        ONLY:  use_seperate_pe_for_dvrp_output
242       
243    USE grid_variables,                                                        &
244        ONLY:  dx, dy
245       
246    USE indices,                                                               &
247        ONLY:  mg_loc_ind, nnx, nny, nnz, nx, ny, nxl_mg, nxr_mg, nyn_mg,      &
248               nys_mg, nzt, nzt_mg
249       
250    USE kinds
[1551]251   
252    USE land_surface_model_mod,                                                &
[1788]253        ONLY:  conserve_water_content, land_surface, nzb_soil,                 &
[1551]254               nzt_soil, root_fraction, soil_moisture, soil_temperature,       &
255               soil_type, soil_type_name, veg_type, veg_type_name, zs
256 
[1320]257    USE model_1d,                                                              &
258        ONLY:  damp_level_ind_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
259       
[1783]260    USE netcdf_interface,                                                      &
261        ONLY:  netcdf_data_format, netcdf_data_format_string, netcdf_deflate
262
[1320]263    USE particle_attributes,                                                   &
264        ONLY:  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, collision_kernel,     &
265               density_ratio, dissipation_classes, dt_min_part, dt_prel,       &
[1359]266               dt_write_particle_data, end_time_prel,                          &
[1320]267               maximum_number_of_tailpoints, maximum_tailpoint_age,            &
268               minimum_tailpoint_distance, number_of_particle_groups,          &
269               particle_advection, particle_advection_start,                   &
270               particles_per_point, pdx, pdy, pdz,  psb, psl, psn, psr, pss,   &
271               pst, radius, radius_classes, random_start_position,             &
[1575]272               seed_follows_topography,                                        &
[1320]273               total_number_of_particles, use_particle_tails,                  &
274               use_sgs_for_particles, total_number_of_tails,                   &
275               vertical_particle_advection, write_particle_statistics
276       
[1]277    USE pegrid
[1484]278
279    USE plant_canopy_model_mod,                                                &
280        ONLY:  alpha_lad, beta_lad, calc_beta_lad_profile, canopy_drag_coeff,  &
281               canopy_mode, cthf, lad, lad_surface, lad_vertical_gradient,     &
282               lad_vertical_gradient_level, lad_vertical_gradient_level_ind,   &
283               lai_beta, leaf_scalar_exch_coeff, leaf_surface_conc, pch_index, &
284               plant_canopy
[1551]285
[1791]286    USE pmc_handle_communicator,                                               &
287        ONLY:  pmc_get_model_info
288
[1764]289    USE pmc_interface,                                                         &
[1791]290        ONLY:  nested_run, nesting_mode
[1764]291
[1551]292    USE radiation_model_mod,                                                   &
[1585]293        ONLY:  albedo, albedo_type, albedo_type_name, constant_albedo,         &
294               day_init, dt_radiation, lambda, lw_radiation, net_radiation,    &
295               radiation, radiation_scheme, sw_radiation, time_utc_init
[1324]296   
297    USE spectrum,                                                              &
298        ONLY:  comp_spectra_level, data_output_sp, plot_spectra_level,         &
299               spectra_direction
[1]300
301    IMPLICIT NONE
302
[1682]303    CHARACTER (LEN=1)  ::  prec                !<
[1320]304   
[1682]305    CHARACTER (LEN=2)  ::  do2d_mode           !<
[1320]306   
[1682]307    CHARACTER (LEN=5)  ::  section_chr         !<
[1320]308   
[1682]309    CHARACTER (LEN=10) ::  coor_chr            !<
310    CHARACTER (LEN=10) ::  host_chr            !<
[1320]311   
[1682]312    CHARACTER (LEN=16) ::  begin_chr           !<
[1320]313   
[1682]314    CHARACTER (LEN=26) ::  ver_rev             !<
[1791]315
316    CHARACTER (LEN=32) ::  cpl_name            !<
[1320]317   
[1682]318    CHARACTER (LEN=40) ::  output_format       !<
[1320]319   
[1682]320    CHARACTER (LEN=70) ::  char1               !<
321    CHARACTER (LEN=70) ::  char2               !<
322    CHARACTER (LEN=70) ::  dopr_chr            !<
323    CHARACTER (LEN=70) ::  do2d_xy             !<
324    CHARACTER (LEN=70) ::  do2d_xz             !<
325    CHARACTER (LEN=70) ::  do2d_yz             !<
326    CHARACTER (LEN=70) ::  do3d_chr            !<
327    CHARACTER (LEN=70) ::  domask_chr          !<
328    CHARACTER (LEN=70) ::  run_classification  !<
[1320]329   
[1682]330    CHARACTER (LEN=85) ::  roben               !<
331    CHARACTER (LEN=85) ::  runten              !<
[1320]332   
[1682]333    CHARACTER (LEN=86) ::  coordinates         !<
334    CHARACTER (LEN=86) ::  gradients           !<
335    CHARACTER (LEN=86) ::  leaf_area_density   !<
336    CHARACTER (LEN=86) ::  roots               !<
337    CHARACTER (LEN=86) ::  slices              !<
338    CHARACTER (LEN=86) ::  temperatures        !<
339    CHARACTER (LEN=86) ::  ugcomponent         !<
340    CHARACTER (LEN=86) ::  vgcomponent         !<
[1]341
[1682]342    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)  !<
[410]343
[1791]344    INTEGER(iwp) ::  av             !<
345    INTEGER(iwp) ::  bh             !<
346    INTEGER(iwp) ::  blx            !<
347    INTEGER(iwp) ::  bly            !<
348    INTEGER(iwp) ::  bxl            !<
349    INTEGER(iwp) ::  bxr            !<
350    INTEGER(iwp) ::  byn            !<
351    INTEGER(iwp) ::  bys            !<
352    INTEGER(iwp) ::  ch             !<
353    INTEGER(iwp) ::  count          !<
354    INTEGER(iwp) ::  cpl_parent_id  !<
355    INTEGER(iwp) ::  cwx            !<
356    INTEGER(iwp) ::  cwy            !<
357    INTEGER(iwp) ::  cxl            !<
358    INTEGER(iwp) ::  cxr            !<
359    INTEGER(iwp) ::  cyn            !<
360    INTEGER(iwp) ::  cys            !<
361    INTEGER(iwp) ::  dim            !<
362    INTEGER(iwp) ::  i              !<
363    INTEGER(iwp) ::  io             !<
364    INTEGER(iwp) ::  j              !<
365    INTEGER(iwp) ::  k              !<
366    INTEGER(iwp) ::  l              !<
367    INTEGER(iwp) ::  ll             !<
368    INTEGER(iwp) ::  mpi_type       !<
369    INTEGER(iwp) ::  my_cpl_id      !<
370    INTEGER(iwp) ::  n              !<
371    INTEGER(iwp) ::  ncpl           !<
372    INTEGER(iwp) ::  npe_total      !<
[1320]373   
[1682]374    REAL(wp) ::  canopy_height                    !< canopy height (in m)
375    REAL(wp) ::  cpuseconds_per_simulated_second  !<
[1791]376    REAL(wp) ::  lower_left_coord_x               !< x-coordinate of nest domain
377    REAL(wp) ::  lower_left_coord_y               !< y-coordinate of nest domain
[1]378
379!
380!-- Open the output file. At the end of the simulation, output is directed
381!-- to unit 19.
382    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
383         .NOT. simulated_time_at_begin /= simulated_time )  THEN
384       io = 15   !  header output on file RUN_CONTROL
385    ELSE
386       io = 19   !  header output on file HEADER
387    ENDIF
388    CALL check_open( io )
389
390!
391!-- At the end of the run, output file (HEADER) will be rewritten with
[1551]392!-- new information
[1]393    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
394
395!
396!-- Determine kind of model run
397    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
[1764]398       run_classification = 'restart run'
[328]399    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[1764]400       run_classification = 'run with cyclic fill of 3D - prerun data'
[147]401    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
[1764]402       run_classification = 'run without 1D - prerun'
[197]403    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
[1764]404       run_classification = 'run with 1D - prerun'
[197]405    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
[1764]406       run_classification = 'run initialized by user'
[1]407    ELSE
[254]408       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
409       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
[1]410    ENDIF
[1764]411    IF ( nested_run )  run_classification = 'nested ' // run_classification
[97]412    IF ( ocean )  THEN
413       run_classification = 'ocean - ' // run_classification
414    ELSE
415       run_classification = 'atmosphere - ' // run_classification
416    ENDIF
[1]417
418!
419!-- Run-identification, date, time, host
420    host_chr = host(1:10)
[75]421    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
[102]422    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
[291]423    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
424#if defined( __mpi2 )
425       mpi_type = 2
426#else
427       mpi_type = 1
428#endif
429       WRITE ( io, 101 )  mpi_type, coupling_mode
430    ENDIF
[1108]431#if defined( __parallel )
[1353]432    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]433       IF ( coupling_start_time > simulated_time_at_begin )  THEN
434          WRITE ( io, 109 )
435       ELSE
436          WRITE ( io, 114 )
437       ENDIF
438    ENDIF
[1108]439#endif
[1429]440    IF ( ensemble_member_nr /= 0 )  THEN
441       WRITE ( io, 512 )  run_date, run_identifier, run_time, runnr,           &
442                       ADJUSTR( host_chr ), ensemble_member_nr
443    ELSE
444       WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr,           &
[102]445                       ADJUSTR( host_chr )
[1429]446    ENDIF
[1]447#if defined( __parallel )
[1482]448    IF ( npex == -1  .AND.  npey == -1 )  THEN
[1]449       char1 = 'calculated'
450    ELSE
451       char1 = 'predefined'
452    ENDIF
453    IF ( threads_per_task == 1 )  THEN
[102]454       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
[1]455    ELSE
[102]456       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
[1]457                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
458    ENDIF
[1111]459    IF ( num_acc_per_node /= 0 )  WRITE ( io, 117 )  num_acc_per_node   
[1]460    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
461           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
462         npex == -1  .AND.  pdims(2) == 1 )                      &
463    THEN
[102]464       WRITE ( io, 106 )
[1]465    ELSEIF ( pdims(2) == 1 )  THEN
[102]466       WRITE ( io, 107 )  'x'
[1]467    ELSEIF ( pdims(1) == 1 )  THEN
[102]468       WRITE ( io, 107 )  'y'
[1]469    ENDIF
[102]470    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
[759]471    IF ( numprocs /= maximum_parallel_io_streams )  THEN
472       WRITE ( io, 108 )  maximum_parallel_io_streams
473    ENDIF
[1111]474#else
475    IF ( num_acc_per_node /= 0 )  WRITE ( io, 120 )  num_acc_per_node
[1]476#endif
[1764]477
478!
479!-- Nesting informations
480    IF ( nested_run )  THEN
[1791]481
482       WRITE ( io, 600 )  TRIM( nesting_mode )
483       CALL pmc_get_model_info( ncpl = ncpl, cpl_id = my_cpl_id )
484
485       DO  n = 1, ncpl
486          CALL pmc_get_model_info( request_for_cpl_id = n, cpl_name = cpl_name,&
487                                   cpl_parent_id = cpl_parent_id,              &
488                                   lower_left_x = lower_left_coord_x,          &
489                                   lower_left_y = lower_left_coord_y,          &
490                                   npe_total = npe_total )
491          IF ( n == my_cpl_id )  THEN
492             char1 = '*'
493          ELSE
494             char1 = ' '
495          ENDIF
496          WRITE ( io, 601 )  TRIM( char1 ), n, cpl_parent_id, npe_total,       &
497                             lower_left_coord_x, lower_left_coord_y,           &
498                             TRIM( cpl_name )
499       ENDDO
[1764]500    ENDIF
[1]501    WRITE ( io, 99 )
502
503!
504!-- Numerical schemes
505    WRITE ( io, 110 )
506    IF ( psolver(1:7) == 'poisfft' )  THEN
507       WRITE ( io, 111 )  TRIM( fft_method )
[1216]508       IF ( transpose_compute_overlap )  WRITE( io, 115 )
[1]509    ELSEIF ( psolver == 'sor' )  THEN
510       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
[1575]511    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
512       WRITE ( io, 135 )  TRIM(psolver), cycle_mg, maximum_grid_level, ngsrb
[1]513       IF ( mg_cycles == -1 )  THEN
514          WRITE ( io, 140 )  residual_limit
515       ELSE
516          WRITE ( io, 141 )  mg_cycles
517       ENDIF
518       IF ( mg_switch_to_pe0_level == 0 )  THEN
519          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
520                             nzt_mg(1)
[197]521       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
[1]522          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
523                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
524                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
525                             nzt_mg(mg_switch_to_pe0_level),    &
526                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
527                             nzt_mg(1)
528       ENDIF
[927]529       IF ( masking_method )  WRITE ( io, 144 )
[1]530    ENDIF
531    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
532    THEN
533       WRITE ( io, 142 )
534    ENDIF
535
536    IF ( momentum_advec == 'pw-scheme' )  THEN
537       WRITE ( io, 113 )
[1299]538    ELSEIF (momentum_advec == 'ws-scheme' )  THEN
[667]539       WRITE ( io, 503 )
[1]540    ENDIF
541    IF ( scalar_advec == 'pw-scheme' )  THEN
542       WRITE ( io, 116 )
[667]543    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
544       WRITE ( io, 504 )
[1557]545    ELSEIF ( scalar_advec == 'ws-scheme-mono' )  THEN
546       WRITE ( io, 513 )
[1]547    ELSE
548       WRITE ( io, 118 )
549    ENDIF
[63]550
551    WRITE ( io, 139 )  TRIM( loop_optimization )
552
[1]553    IF ( galilei_transformation )  THEN
554       IF ( use_ug_for_galilei_tr )  THEN
[868]555          char1 = '0.6 * geostrophic wind'
[1]556       ELSE
557          char1 = 'mean wind in model domain'
558       ENDIF
559       IF ( simulated_time_at_begin == simulated_time )  THEN
560          char2 = 'at the start of the run'
561       ELSE
562          char2 = 'at the end of the run'
563       ENDIF
[1353]564       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ),                        &
565                          advected_distance_x/1000.0_wp,                       &
566                          advected_distance_y/1000.0_wp
[1]567    ENDIF
[1001]568    WRITE ( io, 122 )  timestep_scheme
[87]569    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
[1353]570    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
[108]571       IF ( .NOT. ocean )  THEN
572          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
573               rayleigh_damping_factor
574       ELSE
575          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
576               rayleigh_damping_factor
577       ENDIF
[1]578    ENDIF
[940]579    IF ( neutral )  WRITE ( io, 131 )  pt_surface
[75]580    IF ( humidity )  THEN
[1]581       IF ( .NOT. cloud_physics )  THEN
582          WRITE ( io, 129 )
583       ELSE
584          WRITE ( io, 130 )
585       ENDIF
586    ENDIF
587    IF ( passive_scalar )  WRITE ( io, 134 )
[240]588    IF ( conserve_volume_flow )  THEN
[241]589       WRITE ( io, 150 )  conserve_volume_flow_mode
590       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
591          WRITE ( io, 151 )  u_bulk, v_bulk
592       ENDIF
[240]593    ELSEIF ( dp_external )  THEN
594       IF ( dp_smooth )  THEN
[241]595          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
[240]596       ELSE
[241]597          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
[240]598       ENDIF
599    ENDIF
[1]600    WRITE ( io, 99 )
601
602!
[1551]603!-- Runtime and timestep information
[1]604    WRITE ( io, 200 )
605    IF ( .NOT. dt_fixed )  THEN
606       WRITE ( io, 201 )  dt_max, cfl_factor
607    ELSE
608       WRITE ( io, 202 )  dt
609    ENDIF
610    WRITE ( io, 203 )  simulated_time_at_begin, end_time
611
[1322]612    IF ( time_restart /= 9999999.9_wp  .AND. &
[1]613         simulated_time_at_begin == simulated_time )  THEN
[1322]614       IF ( dt_restart == 9999999.9_wp )  THEN
[1]615          WRITE ( io, 204 )  ' Restart at:       ',time_restart
616       ELSE
617          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
618       ENDIF
619    ENDIF
620
621    IF ( simulated_time_at_begin /= simulated_time )  THEN
622       i = MAX ( log_point_s(10)%counts, 1 )
[1353]623       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0_wp )  THEN
624          cpuseconds_per_simulated_second = 0.0_wp
[1]625       ELSE
626          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
627                                            ( simulated_time -    &
628                                              simulated_time_at_begin )
629       ENDIF
[1322]630       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum,      &
631                          log_point_s(10)%sum / REAL( i, KIND=wp ), &
[1]632                          cpuseconds_per_simulated_second
[1322]633       IF ( time_restart /= 9999999.9_wp  .AND.  time_restart < end_time )  THEN
634          IF ( dt_restart == 9999999.9_wp )  THEN
[1106]635             WRITE ( io, 204 )  ' Next restart at:     ',time_restart
[1]636          ELSE
[1106]637             WRITE ( io, 205 )  ' Next restart at:     ',time_restart, dt_restart
[1]638          ENDIF
639       ENDIF
640    ENDIF
641
[1324]642
[1]643!
[291]644!-- Start time for coupled runs, if independent precursor runs for atmosphere
[1106]645!-- and ocean are used or have been used. In this case, coupling_start_time
646!-- defines the time when the coupling is switched on.
[1353]647    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]648       WRITE ( io, 207 )  coupling_start_time
[291]649    ENDIF
650
651!
[1]652!-- Computational grid
[94]653    IF ( .NOT. ocean )  THEN
654       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
655       IF ( dz_stretch_level_index < nzt+1 )  THEN
656          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
657                             dz_stretch_factor, dz_max
658       ENDIF
659    ELSE
660       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
661       IF ( dz_stretch_level_index > 0 )  THEN
662          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
663                             dz_stretch_factor, dz_max
664       ENDIF
[1]665    ENDIF
666    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
667                       MIN( nnz+2, nzt+2 )
668    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
669
670!
[1365]671!-- Large scale forcing and nudging
672    WRITE ( io, 160 )
673    IF ( large_scale_forcing )  THEN
674       WRITE ( io, 162 )
675       WRITE ( io, 163 )
676
677       IF ( large_scale_subsidence )  THEN
678          IF ( .NOT. use_subsidence_tendencies )  THEN
679             WRITE ( io, 164 )
680          ELSE
681             WRITE ( io, 165 )
682          ENDIF
683       ENDIF
684
685       IF ( bc_pt_b == 'dirichlet' )  THEN
686          WRITE ( io, 180 )
687       ELSEIF ( bc_pt_b == 'neumann' )  THEN
688          WRITE ( io, 181 )
689       ENDIF
690
691       IF ( bc_q_b == 'dirichlet' )  THEN
692          WRITE ( io, 182 )
693       ELSEIF ( bc_q_b == 'neumann' )  THEN
694          WRITE ( io, 183 )
695       ENDIF
696
697       WRITE ( io, 167 )
698       IF ( nudging )  THEN
699          WRITE ( io, 170 )
700       ENDIF
701    ELSE
702       WRITE ( io, 161 )
703       WRITE ( io, 171 )
704    ENDIF
705    IF ( large_scale_subsidence )  THEN
706       WRITE ( io, 168 )
707       WRITE ( io, 169 )
708    ENDIF
709
710!
711!-- Profile for the large scale vertial velocity
712!-- Building output strings, starting with surface value
713    IF ( large_scale_subsidence )  THEN
714       temperatures = '   0.0'
715       gradients = '------'
716       slices = '     0'
717       coordinates = '   0.0'
718       i = 1
719       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
720
721          WRITE (coor_chr,'(E10.2,7X)')  &
722                                w_subs(subs_vertical_gradient_level_i(i))
723          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
724
725          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
726          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
727
728          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
729          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
730
731          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
732          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
733
734          IF ( i == 10 )  THEN
735             EXIT
736          ELSE
737             i = i + 1
738          ENDIF
739
740       ENDDO
741
742 
743       IF ( .NOT. large_scale_forcing )  THEN
744          WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
745                             TRIM( gradients ), TRIM( slices )
746       ENDIF
747
748
749    ENDIF
750
751!-- Profile of the geostrophic wind (component ug)
752!-- Building output strings
753    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
754    gradients = '------'
755    slices = '     0'
756    coordinates = '   0.0'
757    i = 1
758    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
759     
760       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
761       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
762
763       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
764       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
765
766       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
767       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
768
769       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
770       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
771
772       IF ( i == 10 )  THEN
773          EXIT
774       ELSE
775          i = i + 1
776       ENDIF
777
778    ENDDO
779
780    IF ( .NOT. large_scale_forcing )  THEN
781       WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
782                          TRIM( gradients ), TRIM( slices )
783    ENDIF
784
785!-- Profile of the geostrophic wind (component vg)
786!-- Building output strings
787    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
788    gradients = '------'
789    slices = '     0'
790    coordinates = '   0.0'
791    i = 1
792    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
793
794       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
795       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
796
797       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
798       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
799
800       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
801       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
802
803       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
804       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
805
806       IF ( i == 10 )  THEN
807          EXIT
808       ELSE
809          i = i + 1
810       ENDIF
811 
812    ENDDO
813
814    IF ( .NOT. large_scale_forcing )  THEN
815       WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
816                          TRIM( gradients ), TRIM( slices )
817    ENDIF
818
819!
[1]820!-- Topography
821    WRITE ( io, 270 )  topography
822    SELECT CASE ( TRIM( topography ) )
823
824       CASE ( 'flat' )
825          ! no actions necessary
826
827       CASE ( 'single_building' )
828          blx = INT( building_length_x / dx )
829          bly = INT( building_length_y / dy )
[1675]830          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
831          IF ( ABS( zw(bh  ) - building_height ) == &
832               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
[1]833
[1322]834          IF ( building_wall_left == 9999999.9_wp )  THEN
[1]835             building_wall_left = ( nx + 1 - blx ) / 2 * dx
836          ENDIF
[1353]837          bxl = INT ( building_wall_left / dx + 0.5_wp )
[1]838          bxr = bxl + blx
839
[1322]840          IF ( building_wall_south == 9999999.9_wp )  THEN
[1]841             building_wall_south = ( ny + 1 - bly ) / 2 * dy
842          ENDIF
[1353]843          bys = INT ( building_wall_south / dy + 0.5_wp )
[1]844          byn = bys + bly
845
846          WRITE ( io, 271 )  building_length_x, building_length_y, &
847                             building_height, bxl, bxr, bys, byn
848
[240]849       CASE ( 'single_street_canyon' )
[1675]850          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
851          IF ( ABS( zw(ch  ) - canyon_height ) == &
852               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
[1322]853          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[240]854!
855!--          Street canyon in y direction
856             cwx = NINT( canyon_width_x / dx )
[1322]857             IF ( canyon_wall_left == 9999999.9_wp )  THEN
[240]858                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
859             ENDIF
860             cxl = NINT( canyon_wall_left / dx )
861             cxr = cxl + cwx
862             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
863
[1322]864          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[240]865!
866!--          Street canyon in x direction
867             cwy = NINT( canyon_width_y / dy )
[1322]868             IF ( canyon_wall_south == 9999999.9_wp )  THEN
[240]869                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
870             ENDIF
871             cys = NINT( canyon_wall_south / dy )
872             cyn = cys + cwy
873             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
874          ENDIF
875
[1]876    END SELECT
877
[256]878    IF ( TRIM( topography ) /= 'flat' )  THEN
879       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
880          IF ( TRIM( topography ) == 'single_building' .OR.  &
881               TRIM( topography ) == 'single_street_canyon' )  THEN
882             WRITE ( io, 278 )
883          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
884             WRITE ( io, 279 )
885          ENDIF
886       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
887          WRITE ( io, 278 )
888       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
889          WRITE ( io, 279 )
890       ENDIF
891    ENDIF
892
[1299]893    IF ( plant_canopy )  THEN
[1484]894   
895       canopy_height = pch_index * dz
[138]896
[1484]897       WRITE ( io, 280 )  canopy_mode, canopy_height, pch_index,               &
898                          canopy_drag_coeff
[1299]899       IF ( passive_scalar )  THEN
[1484]900          WRITE ( io, 281 )  leaf_scalar_exch_coeff,                           &
901                             leaf_surface_conc
[153]902       ENDIF
[138]903
[1]904!
[153]905!--    Heat flux at the top of vegetation
[1484]906       WRITE ( io, 282 )  cthf
[153]907
908!
[1484]909!--    Leaf area density profile, calculated either from given vertical
910!--    gradients or from beta probability density function.
911       IF (  .NOT.  calc_beta_lad_profile )  THEN
[138]912
[1484]913!--       Building output strings, starting with surface value
[1697]914          WRITE ( leaf_area_density, '(F7.4)' )  lad_surface
[1484]915          gradients = '------'
916          slices = '     0'
917          coordinates = '   0.0'
918          i = 1
919          DO  WHILE ( i < 11  .AND.  lad_vertical_gradient_level_ind(i) /= -9999 )
[138]920
[1484]921             WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
922             leaf_area_density = TRIM( leaf_area_density ) // ' ' // TRIM( coor_chr )
923 
924             WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
925             gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
[138]926
[1484]927             WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
928             slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
[138]929
[1484]930             WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
931             coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
[138]932
[1484]933             i = i + 1
934          ENDDO
[138]935
[1484]936          WRITE ( io, 283 )  TRIM( coordinates ), TRIM( leaf_area_density ),              &
937                             TRIM( gradients ), TRIM( slices )
[138]938
[1484]939       ELSE
940       
[1697]941          WRITE ( leaf_area_density, '(F7.4)' )  lad_surface
[1484]942          coordinates = '   0.0'
943         
944          DO  k = 1, pch_index
945
946             WRITE (coor_chr,'(F7.2)')  lad(k)
947             leaf_area_density = TRIM( leaf_area_density ) // ' ' // TRIM( coor_chr )
948 
949             WRITE (coor_chr,'(F7.1)')  zu(k)
950             coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
951
952          ENDDO       
953
954          WRITE ( io, 284 ) TRIM( coordinates ), TRIM( leaf_area_density ), alpha_lad,    &
955                            beta_lad, lai_beta
956
957       ENDIF 
958
[138]959    ENDIF
960
[1484]961
[1551]962    IF ( land_surface )  THEN
963
964       temperatures = ''
965       gradients    = '' ! use for humidity here
966       coordinates  = '' ! use for height
967       roots        = '' ! use for root fraction
968       slices       = '' ! use for index
969
970       i = 1
971       DO i = nzb_soil, nzt_soil
972          WRITE (coor_chr,'(F10.2,7X)') soil_temperature(i)
973          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
974
975          WRITE (coor_chr,'(F10.2,7X)') soil_moisture(i)
976          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
977
978          WRITE (coor_chr,'(F10.2,7X)')  - zs(i)
979          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
980
981          WRITE (coor_chr,'(F10.2,7X)')  root_fraction(i)
982          roots = TRIM( roots ) // ' '  // TRIM( coor_chr )
983
984          WRITE (coor_chr,'(I10,7X)')  i
985          slices = TRIM( slices ) // ' '  // TRIM( coor_chr )
986
987
988       ENDDO
989
[138]990!
[1551]991!--    Write land surface model header
992       WRITE( io, 419 )
993       IF ( conserve_water_content )  THEN
994          WRITE( io, 440 )
995       ELSE
996          WRITE( io, 441 )
997       ENDIF
998
[1590]999       WRITE( io, 438 ) TRIM( veg_type_name(veg_type) ),                       &
1000                        TRIM (soil_type_name(soil_type) )
[1551]1001       WRITE( io, 439 ) TRIM( coordinates ), TRIM( temperatures ),             &
1002                        TRIM( gradients ), TRIM( roots ), TRIM( slices )
1003
1004
1005    ENDIF
1006
1007    IF ( radiation )  THEN
1008!
[1585]1009!--    Write radiation model header
[1551]1010       WRITE( io, 444 )
1011
1012       IF ( radiation_scheme == "constant" )  THEN
1013          WRITE( io, 445 ) net_radiation
1014       ELSEIF ( radiation_scheme == "clear-sky" )  THEN
1015          WRITE( io, 446 )
[1585]1016       ELSEIF ( radiation_scheme == "rrtmg" )  THEN
1017          WRITE( io, 447 )
1018          IF ( .NOT. lw_radiation )  WRITE( io, 458 )
1019          IF ( .NOT. sw_radiation )  WRITE( io, 459 )
1020       ENDIF
1021
1022       IF ( albedo_type == 0 )  THEN
1023          WRITE( io, 448 ) albedo
[1551]1024       ELSE
[1590]1025          WRITE( io, 456 ) TRIM( albedo_type_name(albedo_type) )
[1551]1026       ENDIF
[1585]1027       IF ( constant_albedo )  THEN
1028          WRITE( io, 457 )
1029       ENDIF
[1551]1030       WRITE( io, 449 ) dt_radiation
1031    ENDIF
1032
1033
1034!
[1]1035!-- Boundary conditions
1036    IF ( ibc_p_b == 0 )  THEN
1037       runten = 'p(0)     = 0      |'
1038    ELSEIF ( ibc_p_b == 1 )  THEN
1039       runten = 'p(0)     = p(1)   |'
1040    ENDIF
1041    IF ( ibc_p_t == 0 )  THEN
1042       roben  = 'p(nzt+1) = 0      |'
1043    ELSE
1044       roben  = 'p(nzt+1) = p(nzt) |'
1045    ENDIF
1046
1047    IF ( ibc_uv_b == 0 )  THEN
1048       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
1049    ELSE
1050       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
1051    ENDIF
[132]1052    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
1053       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
1054    ELSEIF ( ibc_uv_t == 0 )  THEN
[1]1055       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
1056    ELSE
1057       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
1058    ENDIF
1059
1060    IF ( ibc_pt_b == 0 )  THEN
[1551]1061       IF ( land_surface )  THEN
1062          runten = TRIM( runten ) // ' pt(0)     = from soil model'
1063       ELSE
1064          runten = TRIM( runten ) // ' pt(0)     = pt_surface'
1065       ENDIF
[102]1066    ELSEIF ( ibc_pt_b == 1 )  THEN
[1551]1067       runten = TRIM( runten ) // ' pt(0)     = pt(1)'
[102]1068    ELSEIF ( ibc_pt_b == 2 )  THEN
[1551]1069       runten = TRIM( runten ) // ' pt(0)     = from coupled model'
[1]1070    ENDIF
1071    IF ( ibc_pt_t == 0 )  THEN
[19]1072       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
1073    ELSEIF( ibc_pt_t == 1 )  THEN
1074       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
1075    ELSEIF( ibc_pt_t == 2 )  THEN
1076       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
[667]1077
[1]1078    ENDIF
1079
1080    WRITE ( io, 300 )  runten, roben
1081
1082    IF ( .NOT. constant_diffusion )  THEN
1083       IF ( ibc_e_b == 1 )  THEN
1084          runten = 'e(0)     = e(1)'
1085       ELSE
1086          runten = 'e(0)     = e(1) = (u*/0.1)**2'
1087       ENDIF
1088       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
1089
[97]1090       WRITE ( io, 301 )  'e', runten, roben       
[1]1091
1092    ENDIF
1093
[97]1094    IF ( ocean )  THEN
1095       runten = 'sa(0)    = sa(1)'
1096       IF ( ibc_sa_t == 0 )  THEN
1097          roben =  'sa(nzt+1) = sa_surface'
[1]1098       ELSE
[97]1099          roben =  'sa(nzt+1) = sa(nzt)'
[1]1100       ENDIF
[97]1101       WRITE ( io, 301 ) 'sa', runten, roben
1102    ENDIF
[1]1103
[97]1104    IF ( humidity )  THEN
1105       IF ( ibc_q_b == 0 )  THEN
[1551]1106          IF ( land_surface )  THEN
1107             runten = 'q(0)     = from soil model'
1108          ELSE
1109             runten = 'q(0)     = q_surface'
1110          ENDIF
1111
[97]1112       ELSE
1113          runten = 'q(0)     = q(1)'
1114       ENDIF
1115       IF ( ibc_q_t == 0 )  THEN
1116          roben =  'q(nzt)   = q_top'
1117       ELSE
1118          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
1119       ENDIF
1120       WRITE ( io, 301 ) 'q', runten, roben
1121    ENDIF
[1]1122
[97]1123    IF ( passive_scalar )  THEN
1124       IF ( ibc_q_b == 0 )  THEN
1125          runten = 's(0)     = s_surface'
1126       ELSE
1127          runten = 's(0)     = s(1)'
1128       ENDIF
1129       IF ( ibc_q_t == 0 )  THEN
1130          roben =  's(nzt)   = s_top'
1131       ELSE
1132          roben =  's(nzt)   = s(nzt-1) + ds/dz'
1133       ENDIF
1134       WRITE ( io, 301 ) 's', runten, roben
[1]1135    ENDIF
1136
1137    IF ( use_surface_fluxes )  THEN
1138       WRITE ( io, 303 )
1139       IF ( constant_heatflux )  THEN
[1299]1140          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[1241]1141             WRITE ( io, 306 )  shf(0,0)
1142          ELSE
1143             WRITE ( io, 306 )  surface_heatflux
1144          ENDIF
[1]1145          IF ( random_heatflux )  WRITE ( io, 307 )
1146       ENDIF
[75]1147       IF ( humidity  .AND.  constant_waterflux )  THEN
[1299]1148          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[1241]1149             WRITE ( io, 311 ) qsws(0,0)
1150          ELSE
1151             WRITE ( io, 311 ) surface_waterflux
1152          ENDIF
[1]1153       ENDIF
1154       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
1155          WRITE ( io, 313 ) surface_waterflux
1156       ENDIF
1157    ENDIF
1158
[19]1159    IF ( use_top_fluxes )  THEN
1160       WRITE ( io, 304 )
[102]1161       IF ( coupling_mode == 'uncoupled' )  THEN
[151]1162          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
[102]1163          IF ( constant_top_heatflux )  THEN
1164             WRITE ( io, 306 )  top_heatflux
1165          ENDIF
1166       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1167          WRITE ( io, 316 )
[19]1168       ENDIF
[97]1169       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
1170          WRITE ( io, 309 )  top_salinityflux
1171       ENDIF
[75]1172       IF ( humidity  .OR.  passive_scalar )  THEN
[19]1173          WRITE ( io, 315 )
1174       ENDIF
1175    ENDIF
1176
[1691]1177    IF ( constant_flux_layer )  THEN
1178       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length,                     &
1179                          z0h_factor*roughness_length, kappa,                  &
1180                          zeta_min, zeta_max
[1]1181       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
[75]1182       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
[1]1183          WRITE ( io, 312 )
1184       ENDIF
1185       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
1186          WRITE ( io, 314 )
1187       ENDIF
1188    ELSE
1189       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
[1691]1190          WRITE ( io, 310 )  zeta_min, zeta_max
[1]1191       ENDIF
1192    ENDIF
1193
1194    WRITE ( io, 317 )  bc_lr, bc_ns
[707]1195    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1159]1196       WRITE ( io, 318 )  use_cmax, pt_damping_width, pt_damping_factor       
[151]1197       IF ( turbulent_inflow )  THEN
[1560]1198          IF ( .NOT. recycling_yshift ) THEN
1199             WRITE ( io, 319 )  recycling_width, recycling_plane, &
1200                                inflow_damping_height, inflow_damping_width
1201          ELSE
1202             WRITE ( io, 322 )  recycling_width, recycling_plane, &
1203                                inflow_damping_height, inflow_damping_width
1204          END IF
[151]1205       ENDIF
[1]1206    ENDIF
1207
1208!
[1365]1209!-- Initial Profiles
1210    WRITE ( io, 321 )
1211!
1212!-- Initial wind profiles
1213    IF ( u_profile(1) /= 9999999.9_wp )  WRITE ( io, 427 )
1214
1215!
1216!-- Initial temperature profile
1217!-- Building output strings, starting with surface temperature
1218    WRITE ( temperatures, '(F6.2)' )  pt_surface
1219    gradients = '------'
1220    slices = '     0'
1221    coordinates = '   0.0'
1222    i = 1
1223    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1224
1225       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1226       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1227
1228       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1229       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1230
1231       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1232       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1233
1234       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1235       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1236
1237       IF ( i == 10 )  THEN
1238          EXIT
1239       ELSE
1240          i = i + 1
1241       ENDIF
1242
1243    ENDDO
1244
1245    IF ( .NOT. nudging )  THEN
1246       WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1247                          TRIM( gradients ), TRIM( slices )
1248    ELSE
1249       WRITE ( io, 428 ) 
1250    ENDIF
1251
1252!
1253!-- Initial humidity profile
1254!-- Building output strings, starting with surface humidity
1255    IF ( humidity  .OR.  passive_scalar )  THEN
1256       WRITE ( temperatures, '(E8.1)' )  q_surface
1257       gradients = '--------'
1258       slices = '       0'
1259       coordinates = '     0.0'
1260       i = 1
1261       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1262         
1263          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1264          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1265
1266          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1267          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1268         
1269          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1270          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1271         
1272          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1273          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1274
1275          IF ( i == 10 )  THEN
1276             EXIT
1277          ELSE
1278             i = i + 1
1279          ENDIF
1280
1281       ENDDO
1282
1283       IF ( humidity )  THEN
1284          IF ( .NOT. nudging )  THEN
1285             WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1286                                TRIM( gradients ), TRIM( slices )
1287          ENDIF
1288       ELSE
1289          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1290                             TRIM( gradients ), TRIM( slices )
1291       ENDIF
1292    ENDIF
1293
1294!
1295!-- Initial salinity profile
1296!-- Building output strings, starting with surface salinity
1297    IF ( ocean )  THEN
1298       WRITE ( temperatures, '(F6.2)' )  sa_surface
1299       gradients = '------'
1300       slices = '     0'
1301       coordinates = '   0.0'
1302       i = 1
1303       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1304
1305          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1306          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1307
1308          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1309          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1310
1311          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1312          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1313
1314          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1315          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1316
1317          IF ( i == 10 )  THEN
1318             EXIT
1319          ELSE
1320             i = i + 1
1321          ENDIF
1322
1323       ENDDO
1324
1325       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1326                          TRIM( gradients ), TRIM( slices )
1327    ENDIF
1328
1329
1330!
[1]1331!-- Listing of 1D-profiles
[151]1332    WRITE ( io, 325 )  dt_dopr_listing
[1353]1333    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1334       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1335    ENDIF
1336
1337!
1338!-- DATA output
1339    WRITE ( io, 330 )
[1353]1340    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1341       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1342    ENDIF
1343
1344!
1345!-- 1D-profiles
[346]1346    dopr_chr = 'Profile:'
[1]1347    IF ( dopr_n /= 0 )  THEN
1348       WRITE ( io, 331 )
1349
1350       output_format = ''
[1783]1351       output_format = netcdf_data_format_string
1352       IF ( netcdf_deflate == 0 )  THEN
1353          WRITE ( io, 344 )  output_format
1354       ELSE
1355          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1356       ENDIF
[1]1357
1358       DO  i = 1, dopr_n
1359          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
1360          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
1361             WRITE ( io, 332 )  dopr_chr
1362             dopr_chr = '       :'
1363          ENDIF
1364       ENDDO
1365
1366       IF ( dopr_chr /= '' )  THEN
1367          WRITE ( io, 332 )  dopr_chr
1368       ENDIF
1369       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
[1353]1370       IF ( skip_time_dopr /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dopr
[1]1371    ENDIF
1372
1373!
1374!-- 2D-arrays
1375    DO  av = 0, 1
1376
1377       i = 1
1378       do2d_xy = ''
1379       do2d_xz = ''
1380       do2d_yz = ''
1381       DO  WHILE ( do2d(av,i) /= ' ' )
1382
1383          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
1384          do2d_mode = do2d(av,i)(l-1:l)
1385
1386          SELECT CASE ( do2d_mode )
1387             CASE ( 'xy' )
1388                ll = LEN_TRIM( do2d_xy )
1389                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1390             CASE ( 'xz' )
1391                ll = LEN_TRIM( do2d_xz )
1392                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1393             CASE ( 'yz' )
1394                ll = LEN_TRIM( do2d_yz )
1395                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1396          END SELECT
1397
1398          i = i + 1
1399
1400       ENDDO
1401
1402       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
1403              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
[1327]1404              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) ) )  THEN
[1]1405
1406          IF (  av == 0 )  THEN
1407             WRITE ( io, 334 )  ''
1408          ELSE
1409             WRITE ( io, 334 )  '(time-averaged)'
1410          ENDIF
1411
1412          IF ( do2d_at_begin )  THEN
1413             begin_chr = 'and at the start'
1414          ELSE
1415             begin_chr = ''
1416          ENDIF
1417
1418          output_format = ''
[1783]1419          output_format = netcdf_data_format_string
1420          IF ( netcdf_deflate == 0 )  THEN
1421             WRITE ( io, 344 )  output_format
1422          ELSE
1423             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1424          ENDIF
[1]1425
1426          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
1427             i = 1
1428             slices = '/'
1429             coordinates = '/'
1430!
[1551]1431!--          Building strings with index and coordinate information of the
[1]1432!--          slices
1433             DO  WHILE ( section(i,1) /= -9999 )
1434
1435                WRITE (section_chr,'(I5)')  section(i,1)
1436                section_chr = ADJUSTL( section_chr )
1437                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1438
[206]1439                IF ( section(i,1) == -1 )  THEN
[1353]1440                   WRITE (coor_chr,'(F10.1)')  -1.0_wp
[206]1441                ELSE
1442                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
1443                ENDIF
[1]1444                coor_chr = ADJUSTL( coor_chr )
1445                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1446
1447                i = i + 1
1448             ENDDO
1449             IF ( av == 0 )  THEN
1450                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
1451                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
1452                                   TRIM( coordinates )
[1353]1453                IF ( skip_time_do2d_xy /= 0.0_wp )  THEN
[1]1454                   WRITE ( io, 339 )  skip_time_do2d_xy
1455                ENDIF
1456             ELSE
1457                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
1458                                   TRIM( begin_chr ), averaging_interval, &
1459                                   dt_averaging_input, 'k', TRIM( slices ), &
1460                                   TRIM( coordinates )
[1353]1461                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1462                   WRITE ( io, 339 )  skip_time_data_output_av
1463                ENDIF
1464             ENDIF
[1308]1465             IF ( netcdf_data_format > 4 )  THEN
1466                WRITE ( io, 352 )  ntdim_2d_xy(av)
1467             ELSE
1468                WRITE ( io, 353 )
1469             ENDIF
[1]1470          ENDIF
1471
1472          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
1473             i = 1
1474             slices = '/'
1475             coordinates = '/'
1476!
[1551]1477!--          Building strings with index and coordinate information of the
[1]1478!--          slices
1479             DO  WHILE ( section(i,2) /= -9999 )
1480
1481                WRITE (section_chr,'(I5)')  section(i,2)
1482                section_chr = ADJUSTL( section_chr )
1483                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1484
1485                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
1486                coor_chr = ADJUSTL( coor_chr )
1487                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1488
1489                i = i + 1
1490             ENDDO
1491             IF ( av == 0 )  THEN
1492                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
1493                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
1494                                   TRIM( coordinates )
[1353]1495                IF ( skip_time_do2d_xz /= 0.0_wp )  THEN
[1]1496                   WRITE ( io, 339 )  skip_time_do2d_xz
1497                ENDIF
1498             ELSE
1499                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
1500                                   TRIM( begin_chr ), averaging_interval, &
1501                                   dt_averaging_input, 'j', TRIM( slices ), &
1502                                   TRIM( coordinates )
[1353]1503                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1504                   WRITE ( io, 339 )  skip_time_data_output_av
1505                ENDIF
1506             ENDIF
[1308]1507             IF ( netcdf_data_format > 4 )  THEN
1508                WRITE ( io, 352 )  ntdim_2d_xz(av)
1509             ELSE
1510                WRITE ( io, 353 )
1511             ENDIF
[1]1512          ENDIF
1513
1514          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
1515             i = 1
1516             slices = '/'
1517             coordinates = '/'
1518!
[1551]1519!--          Building strings with index and coordinate information of the
[1]1520!--          slices
1521             DO  WHILE ( section(i,3) /= -9999 )
1522
1523                WRITE (section_chr,'(I5)')  section(i,3)
1524                section_chr = ADJUSTL( section_chr )
1525                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1526
1527                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
1528                coor_chr = ADJUSTL( coor_chr )
1529                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1530
1531                i = i + 1
1532             ENDDO
1533             IF ( av == 0 )  THEN
1534                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
1535                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
1536                                   TRIM( coordinates )
[1353]1537                IF ( skip_time_do2d_yz /= 0.0_wp )  THEN
[1]1538                   WRITE ( io, 339 )  skip_time_do2d_yz
1539                ENDIF
1540             ELSE
1541                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
1542                                   TRIM( begin_chr ), averaging_interval, &
1543                                   dt_averaging_input, 'i', TRIM( slices ), &
1544                                   TRIM( coordinates )
[1353]1545                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1546                   WRITE ( io, 339 )  skip_time_data_output_av
1547                ENDIF
1548             ENDIF
[1308]1549             IF ( netcdf_data_format > 4 )  THEN
1550                WRITE ( io, 352 )  ntdim_2d_yz(av)
1551             ELSE
1552                WRITE ( io, 353 )
1553             ENDIF
[1]1554          ENDIF
1555
1556       ENDIF
1557
1558    ENDDO
1559
1560!
1561!-- 3d-arrays
1562    DO  av = 0, 1
1563
1564       i = 1
1565       do3d_chr = ''
1566       DO  WHILE ( do3d(av,i) /= ' ' )
1567
1568          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
1569          i = i + 1
1570
1571       ENDDO
1572
1573       IF ( do3d_chr /= '' )  THEN
1574          IF ( av == 0 )  THEN
1575             WRITE ( io, 336 )  ''
1576          ELSE
1577             WRITE ( io, 336 )  '(time-averaged)'
1578          ENDIF
1579
[1783]1580          output_format = netcdf_data_format_string
1581          IF ( netcdf_deflate == 0 )  THEN
1582             WRITE ( io, 344 )  output_format
1583          ELSE
1584             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1585          ENDIF
[1]1586
1587          IF ( do3d_at_begin )  THEN
1588             begin_chr = 'and at the start'
1589          ELSE
1590             begin_chr = ''
1591          ENDIF
1592          IF ( av == 0 )  THEN
1593             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
1594                                zu(nz_do3d), nz_do3d
1595          ELSE
1596             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
1597                                TRIM( begin_chr ), averaging_interval, &
1598                                dt_averaging_input, zu(nz_do3d), nz_do3d
1599          ENDIF
1600
[1308]1601          IF ( netcdf_data_format > 4 )  THEN
1602             WRITE ( io, 352 )  ntdim_3d(av)
1603          ELSE
1604             WRITE ( io, 353 )
1605          ENDIF
1606
[1]1607          IF ( av == 0 )  THEN
[1353]1608             IF ( skip_time_do3d /= 0.0_wp )  THEN
[1]1609                WRITE ( io, 339 )  skip_time_do3d
1610             ENDIF
1611          ELSE
[1353]1612             IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1613                WRITE ( io, 339 )  skip_time_data_output_av
1614             ENDIF
1615          ENDIF
1616
1617       ENDIF
1618
1619    ENDDO
1620
1621!
[410]1622!-- masked arrays
1623    IF ( masks > 0 )  WRITE ( io, 345 )  &
1624         mask_scale_x, mask_scale_y, mask_scale_z
1625    DO  mid = 1, masks
1626       DO  av = 0, 1
1627
1628          i = 1
1629          domask_chr = ''
1630          DO  WHILE ( domask(mid,av,i) /= ' ' )
1631             domask_chr = TRIM( domask_chr ) // ' ' //  &
1632                          TRIM( domask(mid,av,i) ) // ','
1633             i = i + 1
1634          ENDDO
1635
1636          IF ( domask_chr /= '' )  THEN
1637             IF ( av == 0 )  THEN
1638                WRITE ( io, 346 )  '', mid
1639             ELSE
1640                WRITE ( io, 346 )  ' (time-averaged)', mid
1641             ENDIF
1642
[1783]1643             output_format = netcdf_data_format_string
[1308]1644!--          Parallel output not implemented for mask data, hence
1645!--          output_format must be adjusted.
1646             IF ( netcdf_data_format == 5 ) output_format = 'netCDF4/HDF5'
1647             IF ( netcdf_data_format == 6 ) output_format = 'netCDF4/HDF5 classic'
[1783]1648             IF ( netcdf_deflate == 0 )  THEN
1649                WRITE ( io, 344 )  output_format
1650             ELSE
1651                WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1652             ENDIF
[410]1653
1654             IF ( av == 0 )  THEN
1655                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1656             ELSE
1657                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1658                                   averaging_interval, dt_averaging_input
1659             ENDIF
1660
1661             IF ( av == 0 )  THEN
[1353]1662                IF ( skip_time_domask(mid) /= 0.0_wp )  THEN
[410]1663                   WRITE ( io, 339 )  skip_time_domask(mid)
1664                ENDIF
1665             ELSE
[1353]1666                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[410]1667                   WRITE ( io, 339 )  skip_time_data_output_av
1668                ENDIF
1669             ENDIF
1670!
1671!--          output locations
1672             DO  dim = 1, 3
[1353]1673                IF ( mask(mid,dim,1) >= 0.0_wp )  THEN
[410]1674                   count = 0
[1353]1675                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0_wp )
[410]1676                      count = count + 1
1677                   ENDDO
1678                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1679                                      mask(mid,dim,:count)
[1353]1680                ELSEIF ( mask_loop(mid,dim,1) < 0.0_wp .AND.  &
1681                         mask_loop(mid,dim,2) < 0.0_wp .AND.  &
1682                         mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1683                   WRITE ( io, 350 )  dir(dim), dir(dim)
[1353]1684                ELSEIF ( mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1685                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1686                                      mask_loop(mid,dim,1:2)
1687                ELSE
1688                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1689                                      mask_loop(mid,dim,1:3)
1690                ENDIF
1691             ENDDO
1692          ENDIF
1693
1694       ENDDO
1695    ENDDO
1696
1697!
[1]1698!-- Timeseries
[1322]1699    IF ( dt_dots /= 9999999.9_wp )  THEN
[1]1700       WRITE ( io, 340 )
1701
[1783]1702       output_format = netcdf_data_format_string
1703       IF ( netcdf_deflate == 0 )  THEN
1704          WRITE ( io, 344 )  output_format
1705       ELSE
1706          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1707       ENDIF
[1]1708       WRITE ( io, 341 )  dt_dots
1709    ENDIF
1710
1711#if defined( __dvrp_graphics )
1712!
1713!-- Dvrp-output
[1322]1714    IF ( dt_dvrp /= 9999999.9_wp )  THEN
[1]1715       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1716                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1717       i = 1
1718       l = 0
[336]1719       m = 0
[1]1720       DO WHILE ( mode_dvrp(i) /= ' ' )
1721          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
[130]1722             READ ( mode_dvrp(i), '(10X,I2)' )  j
[1]1723             l = l + 1
1724             IF ( do3d(0,j) /= ' ' )  THEN
[336]1725                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1726                                   isosurface_color(:,l)
[1]1727             ENDIF
1728          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
[130]1729             READ ( mode_dvrp(i), '(6X,I2)' )  j
[336]1730             m = m + 1
1731             IF ( do2d(0,j) /= ' ' )  THEN
1732                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1733                                   slicer_range_limits_dvrp(:,m)
1734             ENDIF
[1]1735          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
[336]1736             WRITE ( io, 363 )  dvrp_psize
1737             IF ( particle_dvrpsize /= 'none' )  THEN
1738                WRITE ( io, 364 )  'size', TRIM( particle_dvrpsize ), &
1739                                   dvrpsize_interval
1740             ENDIF
1741             IF ( particle_color /= 'none' )  THEN
1742                WRITE ( io, 364 )  'color', TRIM( particle_color ), &
1743                                   color_interval
1744             ENDIF
[1]1745          ENDIF
1746          i = i + 1
1747       ENDDO
[237]1748
[336]1749       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1750                          superelevation_y, superelevation, clip_dvrp_l, &
1751                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1752
1753       IF ( TRIM( topography ) /= 'flat' )  THEN
1754          WRITE ( io, 366 )  topography_color
1755          IF ( cluster_size > 1 )  THEN
1756             WRITE ( io, 367 )  cluster_size
1757          ENDIF
[237]1758       ENDIF
1759
[1]1760    ENDIF
1761#endif
1762
1763!
1764!-- Spectra output
[1322]1765    IF ( dt_dosp /= 9999999.9_wp )  THEN
[1]1766       WRITE ( io, 370 )
1767
[1783]1768       output_format = netcdf_data_format_string
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       WRITE ( io, 371 )  dt_dosp
[1353]1775       IF ( skip_time_dosp /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dosp
[1]1776       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
1777                          ( spectra_direction(i), i = 1,10 ),  &
[189]1778                          ( comp_spectra_level(i), i = 1,100 ), &
1779                          ( plot_spectra_level(i), i = 1,100 ), &
[1]1780                          averaging_interval_sp, dt_averaging_input_pr
1781    ENDIF
1782
1783    WRITE ( io, 99 )
1784
1785!
1786!-- Physical quantities
1787    WRITE ( io, 400 )
1788
1789!
1790!-- Geostrophic parameters
[1551]1791    IF ( radiation .AND. radiation_scheme /= 'constant' )  THEN
1792       WRITE ( io, 417 )  lambda
1793    ENDIF
1794    WRITE ( io, 410 )  phi, omega, f, fs
[1]1795
1796!
1797!-- Other quantities
1798    WRITE ( io, 411 )  g
[1551]1799    IF ( radiation .AND. radiation_scheme /= 'constant' )  THEN
1800       WRITE ( io, 418 )  day_init, time_utc_init
1801    ENDIF
1802
[1179]1803    WRITE ( io, 412 )  TRIM( reference_state )
1804    IF ( use_single_reference_value )  THEN
[97]1805       IF ( ocean )  THEN
[1179]1806          WRITE ( io, 413 )  prho_reference
[97]1807       ELSE
[1179]1808          WRITE ( io, 414 )  pt_reference
[97]1809       ENDIF
1810    ENDIF
[1]1811
1812!
1813!-- Cloud physics parameters
[1299]1814    IF ( cloud_physics )  THEN
[57]1815       WRITE ( io, 415 )
1816       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
[1115]1817       IF ( icloud_scheme == 0 )  THEN
[1353]1818          WRITE ( io, 510 ) 1.0E-6_wp * nc_const
[1115]1819          IF ( precipitation )  WRITE ( io, 511 ) c_sedimentation
1820       ENDIF
[1]1821    ENDIF
1822
1823!
[824]1824!-- Cloud physcis parameters / quantities / numerical methods
1825    WRITE ( io, 430 )
1826    IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets)  THEN
1827       WRITE ( io, 431 )
1828    ELSEIF ( humidity  .AND.  cloud_physics )  THEN
1829       WRITE ( io, 432 )
[1496]1830       IF ( cloud_top_radiation )  WRITE ( io, 132 )
[1115]1831       IF ( icloud_scheme == 1 )  THEN
1832          IF ( precipitation )  WRITE ( io, 133 )
1833       ELSEIF ( icloud_scheme == 0 )  THEN
1834          IF ( drizzle )  WRITE ( io, 506 )
1835          IF ( precipitation )  THEN
1836             WRITE ( io, 505 )
1837             IF ( turbulence )  WRITE ( io, 507 )
1838             IF ( ventilation_effect )  WRITE ( io, 508 )
1839             IF ( limiter_sedimentation )  WRITE ( io, 509 )
1840          ENDIF
1841       ENDIF
[824]1842    ELSEIF ( humidity  .AND.  cloud_droplets )  THEN
1843       WRITE ( io, 433 )
1844       IF ( curvature_solution_effects )  WRITE ( io, 434 )
[825]1845       IF ( collision_kernel /= 'none' )  THEN
1846          WRITE ( io, 435 )  TRIM( collision_kernel )
[828]1847          IF ( collision_kernel(6:9) == 'fast' )  THEN
1848             WRITE ( io, 436 )  radius_classes, dissipation_classes
1849          ENDIF
[825]1850       ELSE
[828]1851          WRITE ( io, 437 )
[825]1852       ENDIF
[824]1853    ENDIF
1854
1855!
[1]1856!-- LES / turbulence parameters
1857    WRITE ( io, 450 )
1858
1859!--
1860! ... LES-constants used must still be added here
1861!--
1862    IF ( constant_diffusion )  THEN
1863       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1864                          prandtl_number
1865    ENDIF
1866    IF ( .NOT. constant_diffusion)  THEN
[1353]1867       IF ( e_init > 0.0_wp )  WRITE ( io, 455 )  e_init
1868       IF ( e_min > 0.0_wp )  WRITE ( io, 454 )  e_min
[1]1869       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1870    ENDIF
1871
1872!
1873!-- Special actions during the run
1874    WRITE ( io, 470 )
1875    IF ( create_disturbances )  THEN
1876       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1877                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1878                          zu(disturbance_level_ind_t), disturbance_level_ind_t
[707]1879       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1]1880          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1881       ELSE
1882          WRITE ( io, 473 )  disturbance_energy_limit
1883       ENDIF
1884       WRITE ( io, 474 )  TRIM( random_generator )
1885    ENDIF
[1353]1886    IF ( pt_surface_initial_change /= 0.0_wp )  THEN
[1]1887       WRITE ( io, 475 )  pt_surface_initial_change
1888    ENDIF
[1353]1889    IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1890       WRITE ( io, 476 )  q_surface_initial_change       
1891    ENDIF
[1353]1892    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1893       WRITE ( io, 477 )  q_surface_initial_change       
1894    ENDIF
1895
[60]1896    IF ( particle_advection )  THEN
[1]1897!
[60]1898!--    Particle attributes
1899       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1900                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
[1359]1901                          end_time_prel
[60]1902       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1903       IF ( random_start_position )  WRITE ( io, 481 )
[1575]1904       IF ( seed_follows_topography )  WRITE ( io, 496 )
[60]1905       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1906       WRITE ( io, 495 )  total_number_of_particles
[824]1907       IF ( use_particle_tails  .AND.  maximum_number_of_tailpoints /= 0 )  THEN
[60]1908          WRITE ( io, 483 )  maximum_number_of_tailpoints
1909          IF ( minimum_tailpoint_distance /= 0 )  THEN
1910             WRITE ( io, 484 )  total_number_of_tails,      &
1911                                minimum_tailpoint_distance, &
1912                                maximum_tailpoint_age
1913          ENDIF
[1]1914       ENDIF
[1322]1915       IF ( dt_write_particle_data /= 9999999.9_wp )  THEN
[60]1916          WRITE ( io, 485 )  dt_write_particle_data
[1327]1917          IF ( netcdf_data_format > 1 )  THEN
1918             output_format = 'netcdf (64 bit offset) and binary'
[1]1919          ELSE
[1327]1920             output_format = 'netcdf and binary'
[1]1921          ENDIF
[1783]1922          IF ( netcdf_deflate == 0 )  THEN
1923             WRITE ( io, 344 )  output_format
1924          ELSE
1925             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1926          ENDIF
[1]1927       ENDIF
[1322]1928       IF ( dt_dopts /= 9999999.9_wp )  WRITE ( io, 494 )  dt_dopts
[60]1929       IF ( write_particle_statistics )  WRITE ( io, 486 )
[1]1930
[60]1931       WRITE ( io, 487 )  number_of_particle_groups
[1]1932
[60]1933       DO  i = 1, number_of_particle_groups
[1322]1934          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9_wp )  THEN
[1353]1935             WRITE ( io, 490 )  i, 0.0_wp
[60]1936             WRITE ( io, 492 )
[1]1937          ELSE
[60]1938             WRITE ( io, 490 )  i, radius(i)
[1353]1939             IF ( density_ratio(i) /= 0.0_wp )  THEN
[60]1940                WRITE ( io, 491 )  density_ratio(i)
1941             ELSE
1942                WRITE ( io, 492 )
1943             ENDIF
[1]1944          ENDIF
[60]1945          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1946                             pdx(i), pdy(i), pdz(i)
[336]1947          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
[60]1948       ENDDO
[1]1949
[60]1950    ENDIF
[1]1951
[60]1952
[1]1953!
1954!-- Parameters of 1D-model
1955    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1956       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1957                          mixing_length_1d, dissipation_1d
1958       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1959          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1960       ENDIF
1961    ENDIF
1962
1963!
[1551]1964!-- User-defined information
[1]1965    CALL user_header( io )
1966
1967    WRITE ( io, 99 )
1968
1969!
1970!-- Write buffer contents to disc immediately
[82]1971    CALL local_flush( io )
[1]1972
1973!
1974!-- Here the FORMATs start
1975
1976 99 FORMAT (1X,78('-'))
[1468]1977100 FORMAT (/1X,'******************************',4X,44('-')/        &
1978            1X,'* ',A,' *',4X,A/                               &
1979            1X,'******************************',4X,44('-'))
1980101 FORMAT (35X,'coupled run using MPI-',I1,': ',A/ &
1981            35X,42('-'))
1982102 FORMAT (/' Date:                 ',A8,4X,'Run:       ',A20/      &
1983            ' Time:                 ',A8,4X,'Run-No.:   ',I2.2/     &
[1106]1984            ' Run on host:        ',A10)
[1]1985#if defined( __parallel )
[1468]1986103 FORMAT (' Number of PEs:',10X,I6,4X,'Processor grid (x,y): (',I4,',',I4, &
[1]1987              ')',1X,A)
[1468]1988104 FORMAT (' Number of PEs:',10X,I6,4X,'Tasks:',I4,'   threads per task:',I4/ &
1989              35X,'Processor grid (x,y): (',I4,',',I4,')',1X,A)
1990105 FORMAT (35X,'One additional PE is used to handle'/37X,'the dvrp output!')
1991106 FORMAT (35X,'A 1d-decomposition along x is forced'/ &
1992            35X,'because the job is running on an SMP-cluster')
1993107 FORMAT (35X,'A 1d-decomposition along ',A,' is used')
1994108 FORMAT (35X,'Max. # of parallel I/O streams is ',I5)
1995109 FORMAT (35X,'Precursor run for coupled atmos-ocean run'/ &
1996            35X,42('-'))
1997114 FORMAT (35X,'Coupled atmosphere-ocean run following'/ &
1998            35X,'independent precursor runs'/             &
1999            35X,42('-'))
[1111]2000117 FORMAT (' Accelerator boards / node:  ',I2)
[1]2001#endif
2002110 FORMAT (/' Numerical Schemes:'/ &
2003             ' -----------------'/)
2004111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
2005112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
[1697]2006            '     Iterations (initial/other): ',I3,'/',I3,'  omega =',F6.3)
[1]2007113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
2008                  ' or Upstream')
[1216]2009115 FORMAT ('     FFT and transpositions are overlapping')
[1]2010116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
2011                  ' or Upstream')
2012118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
[1106]2013119 FORMAT (' --> Galilei-Transform applied to horizontal advection:'/ &
2014            '     translation velocity = ',A/ &
[1]2015            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
[1111]2016120 FORMAT (' Accelerator boards: ',8X,I2)
[1]2017122 FORMAT (' --> Time differencing scheme: ',A)
[108]2018123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
[1697]2019            '     maximum damping coefficient:',F6.3, ' 1/s')
[1]2020129 FORMAT (' --> Additional prognostic equation for the specific humidity')
2021130 FORMAT (' --> Additional prognostic equation for the total water content')
[940]2022131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
2023                  F6.2, ' K assumed')
[824]2024132 FORMAT ('     Parameterization of long-wave radiation processes via'/ &
[1]2025            '     effective emissivity scheme')
[824]2026133 FORMAT ('     Precipitation parameterization via Kessler-Scheme')
[1]2027134 FORMAT (' --> Additional prognostic equation for a passive scalar')
[1575]2028135 FORMAT (' --> Solve perturbation pressure via ',A,' method (', &
[1]2029                  A,'-cycle)'/ &
2030            '     number of grid levels:                   ',I2/ &
2031            '     Gauss-Seidel red/black iterations:       ',I2)
2032136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
2033                  I3,')')
2034137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
2035            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
2036                  I3,')'/ &
2037            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
2038                  I3,')')
[63]2039139 FORMAT (' --> Loop optimization method: ',A)
[1]2040140 FORMAT ('     maximum residual allowed:                ',E10.3)
2041141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
2042142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
2043                  'step')
[87]2044143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
2045                  'kinetic energy')
[927]2046144 FORMAT ('     masking method is used')
[1]2047150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
[241]2048                  'conserved'/ &
2049            '     using the ',A,' mode')
2050151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
[306]2051152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
2052           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
2053           /'     starting from dp_level_b =', F8.3, 'm', A /)
[1365]2054160 FORMAT (//' Large scale forcing and nudging:'/ &
2055              ' -------------------------------'/)
2056161 FORMAT (' --> No large scale forcing from external is used (default) ')
2057162 FORMAT (' --> Large scale forcing from external file LSF_DATA is used: ')
2058163 FORMAT ('     - large scale advection tendencies ')
2059164 FORMAT ('     - large scale subsidence velocity w_subs ')
2060165 FORMAT ('     - large scale subsidence tendencies ')
2061167 FORMAT ('     - and geostrophic wind components ug and vg')
2062168 FORMAT (' --> Large-scale vertical motion is used in the ', &
[1299]2063                  'prognostic equation(s) for')
[1365]2064169 FORMAT ('     the scalar(s) only')
2065170 FORMAT (' --> Nudging is used')
2066171 FORMAT (' --> No nudging is used (default) ')
2067180 FORMAT ('     - prescribed surface values for temperature')
[1376]2068181 FORMAT ('     - prescribed surface fluxes for temperature')
2069182 FORMAT ('     - prescribed surface values for humidity')
[1365]2070183 FORMAT ('     - prescribed surface fluxes for humidity')
[1]2071200 FORMAT (//' Run time and time step information:'/ &
2072             ' ----------------------------------'/)
[1106]2073201 FORMAT ( ' Timestep:             variable     maximum value: ',F6.3,' s', &
[1697]2074             '    CFL-factor:',F5.2)
[1106]2075202 FORMAT ( ' Timestep:          dt = ',F6.3,' s'/)
2076203 FORMAT ( ' Start time:          ',F9.3,' s'/ &
2077             ' End time:            ',F9.3,' s')
[1]2078204 FORMAT ( A,F9.3,' s')
2079205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
[1106]2080206 FORMAT (/' Time reached:        ',F9.3,' s'/ &
2081             ' CPU-time used:       ',F9.3,' s     per timestep:               ', &
2082               '  ',F9.3,' s'/                                                    &
[1111]2083             '                                      per second of simulated tim', &
[1]2084               'e: ',F9.3,' s')
[1106]2085207 FORMAT ( ' Coupling start time: ',F9.3,' s')
[1]2086250 FORMAT (//' Computational grid and domain size:'/ &
2087              ' ----------------------------------'// &
2088              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
2089              ' m    dz =    ',F7.3,' m'/ &
2090              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
2091              ' m  z(u) = ',F10.3,' m'/)
2092252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
[1697]2093              ' factor:',F6.3/ &
[1]2094            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
2095254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
2096            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
2097260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
2098             ' degrees')
[1551]2099270 FORMAT (//' Topography information:'/ &
2100              ' ----------------------'// &
[1]2101              1X,'Topography: ',A)
2102271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
2103              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
2104                ' / ',I4)
[240]2105272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
2106              ' direction' / &
2107              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
2108              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
[256]2109278 FORMAT (' Topography grid definition convention:'/ &
2110            ' cell edge (staggered grid points'/  &
2111            ' (u in x-direction, v in y-direction))' /)
2112279 FORMAT (' Topography grid definition convention:'/ &
2113            ' cell center (scalar grid points)' /)
[138]2114280 FORMAT (//' Vegetation canopy (drag) model:'/ &
2115              ' ------------------------------'// &
2116              ' Canopy mode: ', A / &
[1484]2117              ' Canopy height: ',F6.2,'m (',I4,' grid points)' / &
[138]2118              ' Leaf drag coefficient: ',F6.2 /)
[1484]2119281 FORMAT (/ ' Scalar exchange coefficient: ',F6.2 / &
[153]2120              ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /)
2121282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s')
2122283 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
[138]2123              ' Height:              ',A,'  m'/ &
2124              ' Leaf area density:   ',A,'  m**2/m**3'/ &
2125              ' Gradient:            ',A,'  m**2/m**4'/ &
2126              ' Gridpoint:           ',A)
[1484]2127284 FORMAT (//' Characteristic levels of the leaf area density and coefficients:'// &
2128              ' Height:              ',A,'  m'/ &
2129              ' Leaf area density:   ',A,'  m**2/m**3'/ &
2130              ' Coefficient alpha: ',F6.2 / &
2131              ' Coefficient beta: ',F6.2 / &
2132              ' Leaf area index: ',F6.2,'  m**2/m**2' /)
[138]2133               
[1]2134300 FORMAT (//' Boundary conditions:'/ &
2135             ' -------------------'// &
2136             '                     p                    uv             ', &
[1551]2137             '                     pt'// &
[1]2138             ' B. bound.: ',A/ &
2139             ' T. bound.: ',A)
[97]2140301 FORMAT (/'                     ',A// &
[1]2141             ' B. bound.: ',A/ &
2142             ' T. bound.: ',A)
[19]2143303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
2144304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
2145305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
2146               'computational u,v-level:'// &
[1697]2147             '       zp = ',F6.2,' m   z0 =',F7.4,' m   z0h =',F8.5,&
2148             ' m   kappa =',F5.2/ &
2149             '       Rif value range:   ',F8.2,' <= rif <=',F6.2)
[97]2150306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
[1]2151307 FORMAT ('       Heatflux has a random normal distribution')
2152308 FORMAT ('       Predefined surface temperature')
[97]2153309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
[1]2154310 FORMAT (//'    1D-Model:'// &
2155             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
2156311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
2157312 FORMAT ('       Predefined surface humidity')
2158313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
2159314 FORMAT ('       Predefined scalar value at the surface')
[19]2160315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
[102]2161316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
2162                    'atmosphere model')
[1]2163317 FORMAT (//' Lateral boundaries:'/ &
2164            '       left/right:  ',A/    &
2165            '       north/south: ',A)
[1159]2166318 FORMAT (/'       use_cmax: ',L1 / &
2167            '       pt damping layer width = ',F8.2,' m, pt ', &
[1697]2168                    'damping factor =',F7.4)
[151]2169319 FORMAT ('       turbulence recycling at inflow switched on'/ &
2170            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
2171            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
2172320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
[103]2173            '                                          v: ',F9.6,' m**2/s**2')
[1365]2174321 FORMAT (//' Initial profiles:'/ &
2175              ' ----------------')
[1560]2176322 FORMAT ('       turbulence recycling at inflow switched on'/ &
2177            '       y shift of the recycled inflow turbulence switched on'/ &
2178            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
[1592]2179            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m'/)
[151]2180325 FORMAT (//' List output:'/ &
[1]2181             ' -----------'//  &
2182            '    1D-Profiles:'/    &
2183            '       Output every             ',F8.2,' s')
[151]2184326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
[1]2185            '       Averaging input every    ',F8.2,' s')
2186330 FORMAT (//' Data output:'/ &
2187             ' -----------'/)
2188331 FORMAT (/'    1D-Profiles:')
2189332 FORMAT (/'       ',A)
2190333 FORMAT ('       Output every             ',F8.2,' s',/ &
2191            '       Time averaged over       ',F8.2,' s'/ &
2192            '       Averaging input every    ',F8.2,' s')
2193334 FORMAT (/'    2D-Arrays',A,':')
2194335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2195            '       Output every             ',F8.2,' s  ',A/ &
2196            '       Cross sections at ',A1,' = ',A/ &
2197            '       scalar-coordinates:   ',A,' m'/)
2198336 FORMAT (/'    3D-Arrays',A,':')
2199337 FORMAT (/'       Arrays: ',A/ &
2200            '       Output every             ',F8.2,' s  ',A/ &
2201            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
2202339 FORMAT ('       No output during initial ',F8.2,' s')
2203340 FORMAT (/'    Time series:')
2204341 FORMAT ('       Output every             ',F8.2,' s'/)
2205342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2206            '       Output every             ',F8.2,' s  ',A/ &
2207            '       Time averaged over       ',F8.2,' s'/ &
2208            '       Averaging input every    ',F8.2,' s'/ &
2209            '       Cross sections at ',A1,' = ',A/ &
2210            '       scalar-coordinates:   ',A,' m'/)
2211343 FORMAT (/'       Arrays: ',A/ &
2212            '       Output every             ',F8.2,' s  ',A/ &
2213            '       Time averaged over       ',F8.2,' s'/ &
2214            '       Averaging input every    ',F8.2,' s'/ &
2215            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
[292]2216344 FORMAT ('       Output format: ',A/)
[410]2217345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
2218            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
2219            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
2220            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
2221346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
2222347 FORMAT ('       Variables: ',A/ &
2223            '       Output every             ',F8.2,' s')
2224348 FORMAT ('       Variables: ',A/ &
2225            '       Output every             ',F8.2,' s'/ &
2226            '       Time averaged over       ',F8.2,' s'/ &
2227            '       Averaging input every    ',F8.2,' s')
2228349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2229            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
2230            13('       ',8(F8.2,',')/) )
2231350 FORMAT (/'       Output locations in ',A,'-direction: ', &
2232            'all gridpoints along ',A,'-direction (default).' )
2233351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2234            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
2235            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
[1313]2236352 FORMAT  (/'       Number of output time levels allowed: ',I3 /)
2237353 FORMAT  (/'       Number of output time levels allowed: unlimited' /)
[1783]2238354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
[1]2239#if defined( __dvrp_graphics )
2240360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
2241            '       Output every      ',F7.1,' s'/ &
2242            '       Output mode:      ',A/ &
2243            '       Host / User:      ',A,' / ',A/ &
2244            '       Directory:        ',A// &
2245            '       The sequence contains:')
[337]2246361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
2247            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
2248362 FORMAT (/'       Slicer plane ',A/ &
[336]2249            '       Slicer limits: [',F6.2,',',F6.2,']')
[337]2250363 FORMAT (/'       Particles'/ &
[336]2251            '          particle size:  ',F7.2,' m')
2252364 FORMAT ('          particle ',A,' controlled by "',A,'" with interval [', &
2253                       F6.2,',',F6.2,']')
[337]2254365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
[336]2255            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
2256                     ')'/ &
2257            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
2258            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
[337]2259366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
[336]2260367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
[1]2261#endif
2262370 FORMAT ('    Spectra:')
2263371 FORMAT ('       Output every ',F7.1,' s'/)
2264372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
2265            '       Directions: ', 10(A5,',')/                         &
[189]2266            '       height levels  k = ', 20(I3,',')/                  &
2267            '                          ', 20(I3,',')/                  &
2268            '                          ', 20(I3,',')/                  &
2269            '                          ', 20(I3,',')/                  &
2270            '                          ', 19(I3,','),I3,'.'/           &
[1]2271            '       height levels selected for standard plot:'/        &
[189]2272            '                      k = ', 20(I3,',')/                  &
2273            '                          ', 20(I3,',')/                  &
2274            '                          ', 20(I3,',')/                  &
2275            '                          ', 20(I3,',')/                  &
2276            '                          ', 19(I3,','),I3,'.'/           &
[1]2277            '       Time averaged over ', F7.1, ' s,' /                &
2278            '       Profiles for the time averaging are taken every ', &
2279                    F6.1,' s')
2280400 FORMAT (//' Physical quantities:'/ &
2281              ' -------------------'/)
[1551]2282410 FORMAT ('    Geograph. latitude  :   phi    = ',F4.1,' degr'/   &
[1697]2283            '    Angular velocity    :   omega  =',E10.3,' rad/s'/  &
[1551]2284            '    Coriolis parameter  :   f      = ',F9.6,' 1/s'/    &
2285            '                            f*     = ',F9.6,' 1/s')
2286411 FORMAT (/'    Gravity             :   g      = ',F4.1,' m/s**2')
[1179]2287412 FORMAT (/'    Reference state used in buoyancy terms: ',A)
2288413 FORMAT ('       Reference density in buoyancy terms: ',F8.3,' kg/m**3')
2289414 FORMAT ('       Reference temperature in buoyancy terms: ',F8.4,' K')
[1551]2290415 FORMAT (/' Cloud physics parameters:'/ &
2291             ' ------------------------'/)
2292416 FORMAT ('    Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
2293            '    Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
[1697]2294            '    Density of air     :   rho_0 =',F6.3,' kg/m**3'/  &
[1551]2295            '    Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
[1697]2296            '    Vapourization heat :   L_v   =',E9.2,' J/kg')
[1551]2297417 FORMAT ('    Geograph. longitude :   lambda = ',F4.1,' degr')
2298418 FORMAT (/'    Day of the year at model start :   day_init      =     ',I3 &
2299            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
2300419 FORMAT (//' Land surface model information:'/ &
2301              ' ------------------------------'/)
[1]2302420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
2303            '       Height:        ',A,'  m'/ &
2304            '       Temperature:   ',A,'  K'/ &
2305            '       Gradient:      ',A,'  K/100m'/ &
2306            '       Gridpoint:     ',A)
2307421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
2308            '       Height:      ',A,'  m'/ &
2309            '       Humidity:    ',A,'  kg/kg'/ &
2310            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
2311            '       Gridpoint:   ',A)
2312422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
2313            '       Height:                  ',A,'  m'/ &
2314            '       Scalar concentration:    ',A,'  kg/m**3'/ &
2315            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
2316            '       Gridpoint:               ',A)
2317423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
2318            '       Height:      ',A,'  m'/ &
2319            '       ug:          ',A,'  m/s'/ &
2320            '       Gradient:    ',A,'  1/100s'/ &
2321            '       Gridpoint:   ',A)
2322424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
2323            '       Height:      ',A,'  m'/ &
[97]2324            '       vg:          ',A,'  m/s'/ &
[1]2325            '       Gradient:    ',A,'  1/100s'/ &
2326            '       Gridpoint:   ',A)
[97]2327425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
2328            '       Height:     ',A,'  m'/ &
2329            '       Salinity:   ',A,'  psu'/ &
2330            '       Gradient:   ',A,'  psu/100m'/ &
2331            '       Gridpoint:  ',A)
[411]2332426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
2333            '       Height:      ',A,'  m'/ &
2334            '       w_subs:      ',A,'  m/s'/ &
2335            '       Gradient:    ',A,'  (m/s)/100m'/ &
2336            '       Gridpoint:   ',A)
[767]2337427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
2338                  ' profiles')
[1241]2339428 FORMAT (/'    Initial profiles (u, v, pt, q) are taken from file '/ &
2340             '    NUDGING_DATA')
[824]2341430 FORMAT (//' Cloud physics quantities / methods:'/ &
2342              ' ----------------------------------'/)
2343431 FORMAT ('    Humidity is treated as purely passive scalar (no condensati', &
2344                 'on)')
2345432 FORMAT ('    Bulk scheme with liquid water potential temperature and'/ &
2346            '    total water content is used.'/ &
2347            '    Condensation is parameterized via 0% - or 100% scheme.')
2348433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
2349                 'icle model')
2350434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
2351                 ' droplets < 1.0E-6 m')
[825]2352435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
[828]2353436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
2354                    'are used'/ &
2355            '          number of radius classes:       ',I3,'    interval ', &
2356                       '[1.0E-6,2.0E-4] m'/ &
2357            '          number of dissipation classes:   ',I2,'    interval ', &
2358                       '[0,1000] cm**2/s**3')
2359437 FORMAT ('    Droplet collision is switched off')
[1551]2360438 FORMAT (' --> Land surface type  : ',A,/ &
2361            ' --> Soil porosity type : ',A)
2362439 FORMAT (/'    Initial soil temperature and moisture profile:'// &
2363            '       Height:        ',A,'  m'/ &
2364            '       Temperature:   ',A,'  K'/ &
2365            '       Moisture:      ',A,'  m**3/m**3'/ &
2366            '       Root fraction: ',A,'  '/ &
2367            '       Gridpoint:     ',A)
[1788]2368440 FORMAT (' --> Soil bottom is closed (water content is conserved, default)')
2369441 FORMAT (' --> Soil bottom is open (water content is not conserved)')
[1551]2370444 FORMAT (//' Radiation model information:'/                                 &
2371              ' ----------------------------'/)
2372445 FORMAT (' --> Using constant net radiation: net_radiation = ', F6.2, '  W/m**2')
2373446 FORMAT (' --> Simple radiation scheme for clear sky is used (no clouds,',  &
2374                   ' default)')
[1585]2375447 FORMAT (' --> RRTMG scheme is used')
[1697]2376448 FORMAT (/'     User-specific surface albedo: albedo =', F6.3)
[1585]2377449 FORMAT  ('     Timestep: dt_radiation = ', F5.2, '  s')
[1551]2378
[1]2379450 FORMAT (//' LES / Turbulence quantities:'/ &
2380              ' ---------------------------'/)
[824]2381451 FORMAT ('    Diffusion coefficients are constant:'/ &
2382            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
[1697]2383453 FORMAT ('    Mixing length is limited to',F5.2,' * z')
[824]2384454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
2385455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
[1585]2386456 FORMAT (/'    Albedo is set for land surface type: ', A)
2387457 FORMAT (/'    --> Albedo is fixed during the run')
2388458 FORMAT (/'    --> Longwave radiation is disabled')
2389459 FORMAT (/'    --> Shortwave radiation is disabled.')
[1]2390470 FORMAT (//' Actions during the simulation:'/ &
2391              ' -----------------------------'/)
[94]2392471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
[1697]2393            '    Disturbance amplitude           :    ',F5.2, ' m/s'/       &
[94]2394            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
2395            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
[1]2396472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
2397                 ' to i/j =',I4)
2398473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
[1697]2399                 F6.3, ' m**2/s**2')
[1]2400474 FORMAT ('    Random number generator used    : ',A/)
2401475 FORMAT ('    The surface temperature is increased (or decreased, ', &
2402                 'respectively, if'/ &
2403            '    the value is negative) by ',F5.2,' K at the beginning of the',&
2404                 ' 3D-simulation'/)
2405476 FORMAT ('    The surface humidity is increased (or decreased, ',&
2406                 'respectively, if the'/ &
2407            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
2408                 ' the 3D-simulation'/)
2409477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
2410                 'respectively, if the'/ &
2411            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
2412                 ' the 3D-simulation'/)
2413480 FORMAT ('    Particles:'/ &
2414            '    ---------'// &
2415            '       Particle advection is active (switched on at t = ', F7.1, &
2416                    ' s)'/ &
2417            '       Start of new particle generations every  ',F6.1,' s'/ &
2418            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
2419            '                            bottom:     ', A, ' top:         ', A/&
2420            '       Maximum particle age:                 ',F9.1,' s'/ &
[1359]2421            '       Advection stopped at t = ',F9.1,' s'/)
[1]2422481 FORMAT ('       Particles have random start positions'/)
[336]2423482 FORMAT ('          Particles are advected only horizontally'/)
[1]2424483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
2425484 FORMAT ('            Number of tails of the total domain: ',I10/ &
2426            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
2427            '            Maximum age of the end of the tail:  ',F8.2,' s')
2428485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
2429486 FORMAT ('       Particle statistics are written on file'/)
2430487 FORMAT ('       Number of particle groups: ',I2/)
2431488 FORMAT ('       SGS velocity components are used for particle advection'/ &
[1697]2432            '          minimum timestep for advection:', F8.5/)
[1]2433489 FORMAT ('       Number of particles simultaneously released at each ', &
2434                    'point: ', I5/)
2435490 FORMAT ('       Particle group ',I2,':'/ &
2436            '          Particle radius: ',E10.3, 'm')
2437491 FORMAT ('          Particle inertia is activated'/ &
[1697]2438            '             density_ratio (rho_fluid/rho_particle) =',F6.3/)
[1]2439492 FORMAT ('          Particles are advected only passively (no inertia)'/)
2440493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
2441            '                                         y:',F8.1,' - ',F8.1,' m'/&
2442            '                                         z:',F8.1,' - ',F8.1,' m'/&
2443            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
2444                       ' m  dz = ',F8.1,' m'/)
2445494 FORMAT ('       Output of particle time series in NetCDF format every ', &
2446                    F8.2,' s'/)
2447495 FORMAT ('       Number of particles in total domain: ',I10/)
[1575]2448496 FORMAT ('       Initial vertical particle positions are interpreted ', &
2449                    'as relative to the given topography')
[1]2450500 FORMAT (//' 1D-Model parameters:'/                           &
2451              ' -------------------'//                           &
2452            '    Simulation time:                   ',F8.1,' s'/ &
2453            '    Run-controll output every:         ',F8.1,' s'/ &
2454            '    Vertical profile output every:     ',F8.1,' s'/ &
2455            '    Mixing length calculation:         ',A/         &
2456            '    Dissipation calculation:           ',A/)
2457502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
[667]2458503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2459504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
[1115]2460505 FORMAT ('    Precipitation parameterization via Seifert-Beheng-Scheme')
2461506 FORMAT ('    Drizzle parameterization via Stokes law')
2462507 FORMAT ('    Turbulence effects on precipitation process')
2463508 FORMAT ('    Ventilation effects on evaporation of rain drops')
2464509 FORMAT ('    Slope limiter used for sedimentation process')
[1551]2465510 FORMAT ('    Droplet density    :   N_c   = ',F6.1,' 1/cm**3')
2466511 FORMAT ('    Sedimentation Courant number:                  '/&
[1697]2467            '                               C_s   =',F4.1,'        ')
[1429]2468512 FORMAT (/' Date:                 ',A8,6X,'Run:       ',A20/      &
2469            ' Time:                 ',A8,6X,'Run-No.:   ',I2.2/     &
2470            ' Run on host:        ',A10,6X,'En-No.:    ',I2.2)
[1557]2471513 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order ' // & 
2472            '+ monotonic adjustment')
[1791]2473600 FORMAT (/' Nesting informations:'/ &
2474            ' --------------------'/ &
2475            ' Nesting mode:                     ',A// &
2476            ' Nest id  parent  number   lower left coordinates   name'/ &
2477            ' (*=me)     id    of PEs      x (m)     y (m)' )
2478601 FORMAT (2X,A1,1X,I2.2,6X,I2.2,5X,I5,5X,F8.2,2X,F8.2,5X,A)
[1]2479
2480 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.