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

Last change on this file since 3151 was 3083, checked in by gronemeier, 6 years ago

merge with branch rans

  • Property svn:keywords set to Id
File size: 94.9 KB
RevLine 
[1682]1!> @file header.f90
[2000]2!------------------------------------------------------------------------------!
[2696]3! This file is part of the PALM model system.
[1036]4!
[2000]5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
[1036]9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
[2718]17! Copyright 1997-2018 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[254]20! Current revisions:
[1]21! -----------------
[1932]22!
[2701]23!
[1485]24! Former revisions:
25! -----------------
26! $Id: header.f90 3083 2018-06-19 14:03:12Z raasch $
[3083]27! Print RANS-mode constants
28!
29! 3065 2018-06-12 07:03:02Z Giersch
[3065]30! Header output concerning stretching revised
31!
32! 3045 2018-05-28 07:55:41Z Giersch
[3045]33! Error messages revised
34!
35! 2967 2018-04-13 11:22:08Z raasch
[2967]36! bugfix: missing parallel cpp-directives added
37!
38! 2883 2018-03-14 08:29:10Z Giersch
[2883]39! Format of the output of dt_dopr_listing (325) has been changed
40!
41! 2817 2018-02-19 16:32:21Z knoop
[2817]42! Preliminary gust module interface implemented
43!
44! 2776 2018-01-31 10:44:42Z Giersch
[2776]45! Variable synthetic_turbulence_generator has been abbreviated
46!
47! 2746 2018-01-15 12:06:04Z suehring
[2746]48! Move flag plant canopy to modules
49!
50! 2718 2018-01-02 08:49:38Z maronga
[2716]51! Corrected "Former revisions" section
52!
53! 2701 2017-12-15 15:40:50Z suehring
54! Changes from last commit documented
55!
56! 2698 2017-12-14 18:46:24Z suehring
[2701]57! Bugfix in get_topography_top_index
58!
[2716]59! 2696 2017-12-14 17:12:51Z kanani
60! Change in file header (GPL part)
[2696]61! Print information about turbulence closure (TG)
62! Print information about inifor initialization (MS)
63!
64! 2575 2017-10-24 09:57:58Z maronga
[2550]65! Added output for complex terrain simulations
66!
67! 2544 2017-10-13 18:09:32Z maronga
[2544]68! Moved initial day of year and time to inipar.
69!
70! 2339 2017-08-07 13:55:26Z gronemeier
[2339]71! corrected timestamp in header
72!
73! 2338 2017-08-07 12:15:38Z gronemeier
[2338]74! Modularize 1D model
75!
[2339]76! 2320 2017-07-21 12:47:43Z suehring
[2320]77! Modularize large-scale forcing and nudging
78!
79! 2300 2017-06-29 13:31:14Z raasch
[2300]80! host-specific code removed
81!
82! 2299 2017-06-29 10:14:38Z maronga
[2299]83! Modified output for spinups
84!
85! 2298 2017-06-29 09:28:18Z raasch
[2298]86! MPI2 related parts removed
87!
88! 2270 2017-06-09 12:18:47Z maronga
[2270]89! Renamed Prandtl layer to constant flux layer
90!
91! 2259 2017-06-08 09:09:11Z gronemeier
[2259]92! Implemented synthetic turbulence generator
93!
94! 2258 2017-06-08 07:55:13Z suehring
[2258]95! Bugfix, add pre-preprocessor directives to enable non-parrallel mode
96!
97! 2233 2017-05-30 18:08:54Z suehring
[1485]98!
[2233]99! 2232 2017-05-30 17:47:52Z suehring
100! Adjustments to new topography and surface concept
101! Generic tunnel setup added
102!
[2201]103! 2200 2017-04-11 11:37:51Z suehring
104! monotonic_adjustment removed
105!
[2119]106! 2118 2017-01-17 16:38:49Z raasch
107! OpenACC relatec code removed
108!
[2074]109! 2073 2016-11-30 14:34:05Z raasch
110! small bugfix concerning output of scalar profiles
111!
[2051]112! 2050 2016-11-08 15:00:55Z gronemeier
113! Implement turbulent outflow condition
114!
[2038]115! 2037 2016-10-26 11:15:40Z knoop
116! Anelastic approximation implemented
117!
[2001]118! 2000 2016-08-20 18:09:15Z knoop
119! Forced header and separation lines into 80 columns
120!
[1993]121! 1992 2016-08-12 15:14:59Z suehring
122! Adapted for top_scalarflux
123!
[1961]124! 1960 2016-07-12 16:34:24Z suehring
125! Treat humidity and passive scalar separately.
126! Modify misleading information concerning humidity.
127! Bugfix, change unit for humidity flux.
128!
[1958]129! 1957 2016-07-07 10:43:48Z suehring
130! flight module added
131!
[1932]132! 1931 2016-06-10 12:06:59Z suehring
133! Rename multigrid into multigrid_noopt
134!
[1903]135! 1902 2016-05-09 11:18:56Z suehring
136! Write information about masking_method only for multigrid solver
137!
[1851]138! 1849 2016-04-08 11:33:18Z hoffmann
139! Adapted for modularization of microphysics
140!
[1834]141! 1833 2016-04-07 14:23:03Z raasch
142! spectrum renamed spectra_mod, output of spectra related quantities moved to
143! spectra_mod
144!
[1832]145! 1831 2016-04-07 13:15:51Z hoffmann
146! turbulence renamed collision_turbulence,
147! drizzle renamed cloud_water_sedimentation
148!
[1827]149! 1826 2016-04-07 12:01:39Z maronga
150! Moved radiation model header output to the respective module.
151! Moved canopy model header output to the respective module.
152!
[1823]153! 1822 2016-04-07 07:49:42Z hoffmann
154! Tails removed. icloud_scheme replaced by microphysics_*
155!
[1818]156! 1817 2016-04-06 15:44:20Z maronga
157! Moved land_surface_model header output to the respective module.
158!
[1809]159! 1808 2016-04-05 19:44:00Z raasch
160! routine local_flush replaced by FORTRAN statement
161!
[1798]162! 1797 2016-03-21 16:50:28Z raasch
163! output of nesting datatransfer mode
164!
[1792]165! 1791 2016-03-11 10:41:25Z raasch
166! output of nesting informations of all domains
167!
[1789]168! 1788 2016-03-10 11:01:04Z maronga
169! Parameter dewfall removed
170!
[1787]171! 1786 2016-03-08 05:49:27Z raasch
172! cpp-direktives for spectra removed
173!
[1784]174! 1783 2016-03-06 18:36:17Z raasch
175! netcdf module and variable names changed, output of netcdf_deflate
176!
[1765]177! 1764 2016-02-28 12:45:19Z raasch
178! output of nesting informations
179!
[1698]180! 1697 2015-10-28 17:14:10Z raasch
181! small E- and F-FORMAT changes to avoid informative compiler messages about
182! insufficient field width
183!
[1692]184! 1691 2015-10-26 16:17:44Z maronga
185! Renamed prandtl_layer to constant_flux_layer, renames rif_min/rif_max to
186! zeta_min/zeta_max.
187!
[1683]188! 1682 2015-10-07 23:56:08Z knoop
189! Code annotations made doxygen readable
190!
[1676]191! 1675 2015-10-02 08:28:59Z gronemeier
192! Bugfix: Definition of topography grid levels
193!
[1662]194! 1660 2015-09-21 08:15:16Z gronemeier
195! Bugfix: Definition of building/street canyon height if vertical grid stretching
196!         starts below the maximum topography height.
197!
[1591]198! 1590 2015-05-08 13:56:27Z maronga
199! Bugfix: Added TRIM statements for character strings for LSM and radiation code
200!
[1586]201! 1585 2015-04-30 07:05:52Z maronga
202! Further output for radiation model(s).
203!
[1576]204! 1575 2015-03-27 09:56:27Z raasch
205! adjustments for psolver-queries, output of seed_follows_topography
206!
[1561]207! 1560 2015-03-06 10:48:54Z keck
208! output for recycling y shift
209!
[1558]210! 1557 2015-03-05 16:43:04Z suehring
211! output for monotonic limiter
212!
[1552]213! 1551 2015-03-03 14:18:16Z maronga
214! Added informal output for land surface model and radiation model. Removed typo.
215!
[1497]216! 1496 2014-12-02 17:25:50Z maronga
217! Renamed: "radiation -> "cloud_top_radiation"
218!
[1485]219! 1484 2014-10-21 10:53:05Z kanani
[1484]220! Changes due to new module structure of the plant canopy model:
221!   module plant_canopy_model_mod and output for new canopy model parameters
222!   (alpha_lad, beta_lad, lai_beta,...) added,
223!   drag_coefficient, leaf_surface_concentration and scalar_exchange_coefficient
224!   renamed to canopy_drag_coeff, leaf_surface_conc and leaf_scalar_exch_coeff,
225!   learde renamed leaf_area_density.
226! Bugfix: DO-WHILE-loop for lad header information additionally restricted
227! by maximum number of gradient levels (currently 10)
[1483]228!
229! 1482 2014-10-18 12:34:45Z raasch
230! information about calculated or predefined virtual processor topology adjusted
231!
[1469]232! 1468 2014-09-24 14:06:57Z maronga
233! Adapted for use on up to 6-digit processor cores
234!
[1430]235! 1429 2014-07-15 12:53:45Z knoop
236! header exended to provide ensemble_member_nr if specified
237!
[1377]238! 1376 2014-04-26 11:21:22Z boeske
239! Correction of typos
240!
[1366]241! 1365 2014-04-22 15:03:56Z boeske
242! New section 'Large scale forcing and nudging':
243! output of large scale forcing and nudging information,
244! new section for initial profiles created
245!
[1360]246! 1359 2014-04-11 17:15:14Z hoffmann
247! dt_sort_particles removed
248!
[1354]249! 1353 2014-04-08 15:21:23Z heinze
250! REAL constants provided with KIND-attribute
251!
[1329]252! 1327 2014-03-21 11:00:16Z raasch
253! parts concerning iso2d and avs output removed,
254! -netcdf output queries
255!
[1325]256! 1324 2014-03-21 09:13:16Z suehring
257! Bugfix: module spectrum added
258!
[1323]259! 1322 2014-03-20 16:38:49Z raasch
260! REAL functions provided with KIND-attribute,
261! some REAL constants defined as wp-kind
262!
[1321]263! 1320 2014-03-20 08:40:49Z raasch
[1320]264! ONLY-attribute added to USE-statements,
265! kind-parameters added to all INTEGER and REAL declaration statements,
266! kinds are defined in new module kinds,
267! revision history before 2012 removed,
268! comment fields (!:) to be used for variable explanations added to
269! all variable declaration statements
[1321]270!
[1309]271! 1308 2014-03-13 14:58:42Z fricke
272! output of the fixed number of output time levels
273! output_format adjusted for masked data if netcdf_data_format > 5
274!
[1300]275! 1299 2014-03-06 13:15:21Z heinze
276! output for using large_scale subsidence in combination
277! with large_scale_forcing
278! reformatting, more detailed explanations
279!
[1242]280! 1241 2013-10-30 11:36:58Z heinze
281! output for nudging + large scale forcing from external file
282!
[1217]283! 1216 2013-08-26 09:31:42Z raasch
284! output for transpose_compute_overlap
285!
[1213]286! 1212 2013-08-15 08:46:27Z raasch
287! output for poisfft_hybrid removed
288!
[1182]289! 1179 2013-06-14 05:57:58Z raasch
290! output of reference_state, use_reference renamed use_single_reference_value
291!
[1160]292! 1159 2013-05-21 11:58:22Z fricke
293! +use_cmax
294!
[1116]295! 1115 2013-03-26 18:16:16Z hoffmann
296! descriptions for Seifert-Beheng-cloud-physics-scheme added
297!
[1112]298! 1111 2013-03-08 23:54:10Z raasch
299! output of accelerator board information
300! ibc_p_b = 2 removed
301!
[1109]302! 1108 2013-03-05 07:03:32Z raasch
303! bugfix for r1106
304!
[1107]305! 1106 2013-03-04 05:31:38Z raasch
306! some format changes for coupled runs
307!
[1093]308! 1092 2013-02-02 11:24:22Z raasch
309! unused variables removed
310!
[1037]311! 1036 2012-10-22 13:43:42Z raasch
312! code put under GPL (PALM 3.9)
313!
[1035]314! 1031 2012-10-19 14:35:30Z raasch
315! output of netCDF data format modified
316!
[1017]317! 1015 2012-09-27 09:23:24Z raasch
[1365]318! output of Adjustment of mixing length to the Prandtl mixing length at first
[1017]319! grid point above ground removed
320!
[1004]321! 1003 2012-09-14 14:35:53Z raasch
322! output of information about equal/unequal subdomain size removed
323!
[1002]324! 1001 2012-09-13 14:08:46Z raasch
325! all actions concerning leapfrog- and upstream-spline-scheme removed
326!
[979]327! 978 2012-08-09 08:28:32Z fricke
328! -km_damp_max, outflow_damping_width
329! +pt_damping_factor, pt_damping_width
330! +z0h
331!
[965]332! 964 2012-07-26 09:14:24Z raasch
333! output of profil-related quantities removed
334!
[941]335! 940 2012-07-09 14:31:00Z raasch
336! Output in case of simulations for pure neutral stratification (no pt-equation
337! solved)
338!
[928]339! 927 2012-06-06 19:15:04Z raasch
340! output of masking_method for mg-solver
341!
[869]342! 868 2012-03-28 12:21:07Z raasch
343! translation velocity in Galilean transformation changed to 0.6 * ug
344!
[834]345! 833 2012-02-22 08:55:55Z maronga
346! Adjusted format for leaf area density
347!
[829]348! 828 2012-02-21 12:00:36Z raasch
349! output of dissipation_classes + radius_classes
350!
[826]351! 825 2012-02-19 03:03:44Z raasch
352! Output of cloud physics parameters/quantities complemented and restructured
353!
[1]354! Revision 1.1  1997/08/11 06:17:20  raasch
355! Initial revision
356!
357!
358! Description:
359! ------------
[1764]360!> Writing a header with all important information about the current run.
[1682]361!> This subroutine is called three times, two times at the beginning
362!> (writing information on files RUN_CONTROL and HEADER) and one time at the
363!> end of the run, then writing additional information about CPU-usage on file
364!> header.
[411]365!-----------------------------------------------------------------------------!
[1682]366 SUBROUTINE header
367 
[1]368
[1320]369    USE arrays_3d,                                                             &
[2232]370        ONLY:  pt_init, q_init, s_init, sa_init, ug, vg, w_subs, zu, zw
[1320]371       
[1]372    USE control_parameters
[1320]373       
374    USE cloud_parameters,                                                      &
[1849]375        ONLY:  cp, l_v, r_d
376
[1320]377    USE cpulog,                                                                &
378        ONLY:  log_point_s
[2544]379           
380    USE date_and_time_mod,                                                     &
381        ONLY:  day_of_year_init, time_utc_init
382
[1320]383    USE dvrp_variables,                                                        &
384        ONLY:  use_seperate_pe_for_dvrp_output
385       
[1957]386    USE flight_mod,                                                            &
387        ONLY:  flight_header
388       
[1320]389    USE grid_variables,                                                        &
390        ONLY:  dx, dy
[2817]391
392    USE gust_mod,                                                              &
393        ONLY: gust_header, gust_module_enabled
[1320]394       
395    USE indices,                                                               &
396        ONLY:  mg_loc_ind, nnx, nny, nnz, nx, ny, nxl_mg, nxr_mg, nyn_mg,      &
397               nys_mg, nzt, nzt_mg
398       
399    USE kinds
[1817]400 
[1551]401    USE land_surface_model_mod,                                                &
[2232]402        ONLY: lsm_header
[1849]403
[2320]404    USE lsf_nudging_mod,                                                       &
405        ONLY:  lsf_nudging_header
406
[1849]407    USE microphysics_mod,                                                      &
408        ONLY:  cloud_water_sedimentation, collision_turbulence,                &
409               c_sedimentation, limiter_sedimentation, nc_const,               &
410               ventilation_effect
411
[2338]412    USE model_1d_mod,                                                          &
[1320]413        ONLY:  damp_level_ind_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
[2338]414
[1783]415    USE netcdf_interface,                                                      &
416        ONLY:  netcdf_data_format, netcdf_data_format_string, netcdf_deflate
417
[1320]418    USE particle_attributes,                                                   &
419        ONLY:  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, collision_kernel,     &
[1831]420               curvature_solution_effects,                                     &
[1320]421               density_ratio, dissipation_classes, dt_min_part, dt_prel,       &
[1359]422               dt_write_particle_data, end_time_prel,                          &
[1822]423               number_of_particle_groups, particle_advection,                  &
424               particle_advection_start,                                       &
[1320]425               particles_per_point, pdx, pdy, pdz,  psb, psl, psn, psr, pss,   &
426               pst, radius, radius_classes, random_start_position,             &
[1575]427               seed_follows_topography,                                        &
[1822]428               total_number_of_particles, use_sgs_for_particles,               &
[1320]429               vertical_particle_advection, write_particle_statistics
430       
[1]431    USE pegrid
[1484]432
433    USE plant_canopy_model_mod,                                                &
[2746]434        ONLY:  pcm_header
[2259]435
[2967]436#if defined( __parallel )
[1791]437    USE pmc_handle_communicator,                                               &
438        ONLY:  pmc_get_model_info
[2967]439#endif
[1791]440
[1764]441    USE pmc_interface,                                                         &
[1797]442        ONLY:  nested_run, nesting_datatransfer_mode, nesting_mode
[1764]443
[1551]444    USE radiation_model_mod,                                                   &
[1826]445        ONLY:  radiation, radiation_header
[1324]446   
[1833]447    USE spectra_mod,                                                           &
448        ONLY:  calculate_spectra, spectra_header
[1]449
[2232]450    USE surface_mod,                                                           &
[2698]451        ONLY:  surf_def_h, get_topography_top_index_ji
[2232]452
[2259]453    USE synthetic_turbulence_generator_mod,                                    &
454        ONLY:  stg_header
455
[3083]456    USE turbulence_closure_mod,                                                &
457        ONLY:  rans_const_c, rans_const_sigma
458
[1]459    IMPLICIT NONE
460
[1682]461    CHARACTER (LEN=1)  ::  prec                !<
[1320]462   
[1682]463    CHARACTER (LEN=2)  ::  do2d_mode           !<
[1320]464   
[1682]465    CHARACTER (LEN=5)  ::  section_chr         !<
[1320]466   
[1682]467    CHARACTER (LEN=10) ::  coor_chr            !<
468    CHARACTER (LEN=10) ::  host_chr            !<
[1320]469   
[1682]470    CHARACTER (LEN=16) ::  begin_chr           !<
[1320]471   
[1682]472    CHARACTER (LEN=26) ::  ver_rev             !<
[1791]473
474    CHARACTER (LEN=32) ::  cpl_name            !<
[1320]475   
[1682]476    CHARACTER (LEN=40) ::  output_format       !<
[1320]477   
[1682]478    CHARACTER (LEN=70) ::  char1               !<
479    CHARACTER (LEN=70) ::  char2               !<
480    CHARACTER (LEN=70) ::  dopr_chr            !<
481    CHARACTER (LEN=70) ::  do2d_xy             !<
482    CHARACTER (LEN=70) ::  do2d_xz             !<
483    CHARACTER (LEN=70) ::  do2d_yz             !<
484    CHARACTER (LEN=70) ::  do3d_chr            !<
485    CHARACTER (LEN=70) ::  domask_chr          !<
486    CHARACTER (LEN=70) ::  run_classification  !<
[1320]487   
[1826]488    CHARACTER (LEN=85) ::  r_upper             !<
489    CHARACTER (LEN=85) ::  r_lower             !<
[1320]490   
[1682]491    CHARACTER (LEN=86) ::  coordinates         !<
492    CHARACTER (LEN=86) ::  gradients           !<
493    CHARACTER (LEN=86) ::  slices              !<
494    CHARACTER (LEN=86) ::  temperatures        !<
495    CHARACTER (LEN=86) ::  ugcomponent         !<
496    CHARACTER (LEN=86) ::  vgcomponent         !<
[1]497
[1682]498    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)  !<
[410]499
[1791]500    INTEGER(iwp) ::  av             !<
501    INTEGER(iwp) ::  bh             !<
502    INTEGER(iwp) ::  blx            !<
503    INTEGER(iwp) ::  bly            !<
504    INTEGER(iwp) ::  bxl            !<
505    INTEGER(iwp) ::  bxr            !<
506    INTEGER(iwp) ::  byn            !<
507    INTEGER(iwp) ::  bys            !<
508    INTEGER(iwp) ::  ch             !<
509    INTEGER(iwp) ::  count          !<
510    INTEGER(iwp) ::  cpl_parent_id  !<
511    INTEGER(iwp) ::  cwx            !<
512    INTEGER(iwp) ::  cwy            !<
513    INTEGER(iwp) ::  cxl            !<
514    INTEGER(iwp) ::  cxr            !<
515    INTEGER(iwp) ::  cyn            !<
516    INTEGER(iwp) ::  cys            !<
517    INTEGER(iwp) ::  dim            !<
518    INTEGER(iwp) ::  i              !<
519    INTEGER(iwp) ::  io             !<
520    INTEGER(iwp) ::  j              !<
521    INTEGER(iwp) ::  k              !<
522    INTEGER(iwp) ::  l              !<
523    INTEGER(iwp) ::  ll             !<
524    INTEGER(iwp) ::  my_cpl_id      !<
525    INTEGER(iwp) ::  n              !<
526    INTEGER(iwp) ::  ncpl           !<
527    INTEGER(iwp) ::  npe_total      !<
[1320]528   
[1826]529
[1682]530    REAL(wp) ::  cpuseconds_per_simulated_second  !<
[1791]531    REAL(wp) ::  lower_left_coord_x               !< x-coordinate of nest domain
532    REAL(wp) ::  lower_left_coord_y               !< y-coordinate of nest domain
[1]533
534!
535!-- Open the output file. At the end of the simulation, output is directed
536!-- to unit 19.
537    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
538         .NOT. simulated_time_at_begin /= simulated_time )  THEN
539       io = 15   !  header output on file RUN_CONTROL
540    ELSE
541       io = 19   !  header output on file HEADER
542    ENDIF
543    CALL check_open( io )
544
545!
546!-- At the end of the run, output file (HEADER) will be rewritten with
[1551]547!-- new information
[1]548    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
549
550!
551!-- Determine kind of model run
552    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
[1764]553       run_classification = 'restart run'
[328]554    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[1764]555       run_classification = 'run with cyclic fill of 3D - prerun data'
[147]556    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
[1764]557       run_classification = 'run without 1D - prerun'
[197]558    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
[1764]559       run_classification = 'run with 1D - prerun'
[2696]560    ELSEIF ( INDEX( initializing_actions, 'inifor' ) /= 0 )  THEN
561       run_classification = 'run initialized with COSMO data'
[197]562    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
[1764]563       run_classification = 'run initialized by user'
[3045]564    ELSEIF ( INDEX( initializing_actions, 'initialize_vortex' ) /=0 )  THEN
565       run_classification = 'run additionally initialized by a Rankine-vortex'
566    ELSEIF ( INDEX( initializing_actions, 'initialize_ptanom' ) /=0 )  THEN
567       run_classification = 'run additionally initialized by temperature anomaly'
[1]568    ELSE
[254]569       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
570       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
[1]571    ENDIF
[1764]572    IF ( nested_run )  run_classification = 'nested ' // run_classification
[97]573    IF ( ocean )  THEN
574       run_classification = 'ocean - ' // run_classification
575    ELSE
576       run_classification = 'atmosphere - ' // run_classification
577    ENDIF
[1]578
579!
580!-- Run-identification, date, time, host
581    host_chr = host(1:10)
[75]582    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
[102]583    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
[291]584    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
[2298]585       WRITE ( io, 101 )  coupling_mode
[291]586    ENDIF
[1108]587#if defined( __parallel )
[2299]588    IF ( coupling_start_time /= 0.0_wp  .AND. .NOT. spinup )  THEN
[1106]589       IF ( coupling_start_time > simulated_time_at_begin )  THEN
590          WRITE ( io, 109 )
591       ELSE
592          WRITE ( io, 114 )
593       ENDIF
594    ENDIF
[1108]595#endif
[1429]596    IF ( ensemble_member_nr /= 0 )  THEN
597       WRITE ( io, 512 )  run_date, run_identifier, run_time, runnr,           &
598                       ADJUSTR( host_chr ), ensemble_member_nr
599    ELSE
600       WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr,           &
[102]601                       ADJUSTR( host_chr )
[1429]602    ENDIF
[1]603#if defined( __parallel )
[1482]604    IF ( npex == -1  .AND.  npey == -1 )  THEN
[1]605       char1 = 'calculated'
606    ELSE
607       char1 = 'predefined'
608    ENDIF
609    IF ( threads_per_task == 1 )  THEN
[102]610       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
[1]611    ELSE
[102]612       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
[1]613                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
614    ENDIF
[2300]615
616    IF ( pdims(2) == 1 )  THEN
[102]617       WRITE ( io, 107 )  'x'
[1]618    ELSEIF ( pdims(1) == 1 )  THEN
[102]619       WRITE ( io, 107 )  'y'
[1]620    ENDIF
[102]621    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
[759]622    IF ( numprocs /= maximum_parallel_io_streams )  THEN
623       WRITE ( io, 108 )  maximum_parallel_io_streams
624    ENDIF
[1]625#endif
[1764]626
627!
628!-- Nesting informations
629    IF ( nested_run )  THEN
[1791]630
[2967]631#if defined( __parallel )
[1797]632       WRITE ( io, 600 )  TRIM( nesting_mode ),                                &
633                          TRIM( nesting_datatransfer_mode )
[1791]634       CALL pmc_get_model_info( ncpl = ncpl, cpl_id = my_cpl_id )
635
636       DO  n = 1, ncpl
637          CALL pmc_get_model_info( request_for_cpl_id = n, cpl_name = cpl_name,&
638                                   cpl_parent_id = cpl_parent_id,              &
639                                   lower_left_x = lower_left_coord_x,          &
640                                   lower_left_y = lower_left_coord_y,          &
641                                   npe_total = npe_total )
642          IF ( n == my_cpl_id )  THEN
643             char1 = '*'
644          ELSE
645             char1 = ' '
646          ENDIF
647          WRITE ( io, 601 )  TRIM( char1 ), n, cpl_parent_id, npe_total,       &
648                             lower_left_coord_x, lower_left_coord_y,           &
649                             TRIM( cpl_name )
650       ENDDO
[2967]651#endif
652
[1764]653    ENDIF
[1]654    WRITE ( io, 99 )
655
656!
657!-- Numerical schemes
658    WRITE ( io, 110 )
[2696]659    IF ( rans_mode )  THEN
660       WRITE ( io, 124 )  TRIM( turbulence_closure ), 'RANS'
661    ELSE
662       WRITE ( io, 124 )  TRIM( turbulence_closure ), 'LES'
663    ENDIF
[2037]664    WRITE ( io, 121 )  TRIM( approximation )
[1]665    IF ( psolver(1:7) == 'poisfft' )  THEN
666       WRITE ( io, 111 )  TRIM( fft_method )
[1216]667       IF ( transpose_compute_overlap )  WRITE( io, 115 )
[1]668    ELSEIF ( psolver == 'sor' )  THEN
669       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
[1575]670    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
671       WRITE ( io, 135 )  TRIM(psolver), cycle_mg, maximum_grid_level, ngsrb
[1]672       IF ( mg_cycles == -1 )  THEN
673          WRITE ( io, 140 )  residual_limit
674       ELSE
675          WRITE ( io, 141 )  mg_cycles
676       ENDIF
677       IF ( mg_switch_to_pe0_level == 0 )  THEN
678          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
679                             nzt_mg(1)
[197]680       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
[1]681          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
682                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
683                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
684                             nzt_mg(mg_switch_to_pe0_level),    &
685                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
686                             nzt_mg(1)
687       ENDIF
[1931]688       IF ( psolver == 'multigrid_noopt' .AND. masking_method )  WRITE ( io, 144 )
[1]689    ENDIF
690    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
691    THEN
692       WRITE ( io, 142 )
693    ENDIF
694
695    IF ( momentum_advec == 'pw-scheme' )  THEN
696       WRITE ( io, 113 )
[1299]697    ELSEIF (momentum_advec == 'ws-scheme' )  THEN
[667]698       WRITE ( io, 503 )
[1]699    ENDIF
700    IF ( scalar_advec == 'pw-scheme' )  THEN
701       WRITE ( io, 116 )
[667]702    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
703       WRITE ( io, 504 )
[1]704    ELSE
705       WRITE ( io, 118 )
706    ENDIF
[63]707
708    WRITE ( io, 139 )  TRIM( loop_optimization )
709
[1]710    IF ( galilei_transformation )  THEN
711       IF ( use_ug_for_galilei_tr )  THEN
[868]712          char1 = '0.6 * geostrophic wind'
[1]713       ELSE
714          char1 = 'mean wind in model domain'
715       ENDIF
716       IF ( simulated_time_at_begin == simulated_time )  THEN
717          char2 = 'at the start of the run'
718       ELSE
719          char2 = 'at the end of the run'
720       ENDIF
[1353]721       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ),                        &
722                          advected_distance_x/1000.0_wp,                       &
723                          advected_distance_y/1000.0_wp
[1]724    ENDIF
[1001]725    WRITE ( io, 122 )  timestep_scheme
[87]726    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
[1353]727    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
[108]728       IF ( .NOT. ocean )  THEN
729          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
730               rayleigh_damping_factor
731       ELSE
732          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
733               rayleigh_damping_factor
734       ENDIF
[1]735    ENDIF
[940]736    IF ( neutral )  WRITE ( io, 131 )  pt_surface
[75]737    IF ( humidity )  THEN
[1]738       IF ( .NOT. cloud_physics )  THEN
739          WRITE ( io, 129 )
740       ELSE
741          WRITE ( io, 130 )
742       ENDIF
743    ENDIF
744    IF ( passive_scalar )  WRITE ( io, 134 )
[240]745    IF ( conserve_volume_flow )  THEN
[241]746       WRITE ( io, 150 )  conserve_volume_flow_mode
747       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
748          WRITE ( io, 151 )  u_bulk, v_bulk
749       ENDIF
[240]750    ELSEIF ( dp_external )  THEN
751       IF ( dp_smooth )  THEN
[241]752          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
[240]753       ELSE
[241]754          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
[240]755       ENDIF
756    ENDIF
[1]757    WRITE ( io, 99 )
758
759!
[1551]760!-- Runtime and timestep information
[1]761    WRITE ( io, 200 )
762    IF ( .NOT. dt_fixed )  THEN
763       WRITE ( io, 201 )  dt_max, cfl_factor
764    ELSE
765       WRITE ( io, 202 )  dt
766    ENDIF
767    WRITE ( io, 203 )  simulated_time_at_begin, end_time
768
[1322]769    IF ( time_restart /= 9999999.9_wp  .AND. &
[1]770         simulated_time_at_begin == simulated_time )  THEN
[1322]771       IF ( dt_restart == 9999999.9_wp )  THEN
[1]772          WRITE ( io, 204 )  ' Restart at:       ',time_restart
773       ELSE
774          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
775       ENDIF
776    ENDIF
777
778    IF ( simulated_time_at_begin /= simulated_time )  THEN
779       i = MAX ( log_point_s(10)%counts, 1 )
[1353]780       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0_wp )  THEN
781          cpuseconds_per_simulated_second = 0.0_wp
[1]782       ELSE
783          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
784                                            ( simulated_time -    &
785                                              simulated_time_at_begin )
786       ENDIF
[1322]787       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum,      &
788                          log_point_s(10)%sum / REAL( i, KIND=wp ), &
[1]789                          cpuseconds_per_simulated_second
[1322]790       IF ( time_restart /= 9999999.9_wp  .AND.  time_restart < end_time )  THEN
791          IF ( dt_restart == 9999999.9_wp )  THEN
[1106]792             WRITE ( io, 204 )  ' Next restart at:     ',time_restart
[1]793          ELSE
[1106]794             WRITE ( io, 205 )  ' Next restart at:     ',time_restart, dt_restart
[1]795          ENDIF
796       ENDIF
797    ENDIF
798
[1324]799
[1]800!
[291]801!-- Start time for coupled runs, if independent precursor runs for atmosphere
[1106]802!-- and ocean are used or have been used. In this case, coupling_start_time
803!-- defines the time when the coupling is switched on.
[1353]804    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]805       WRITE ( io, 207 )  coupling_start_time
[291]806    ENDIF
807
808!
[1]809!-- Computational grid
[94]810    IF ( .NOT. ocean )  THEN
[3065]811       WRITE ( io, 250 )  dx, dy
812       
813       DO i = 1, number_stretch_level_start+1
814          WRITE ( io, 253 )  i, dz(i)
815       ENDDO
816       
817       WRITE( io, 251 ) (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
818       
819       IF ( ANY( dz_stretch_level_start_index < nzt+1 ) )  THEN
820          WRITE( io, '(A)', advance='no') ' Vertical stretching starts at height:'
821          DO i = 1, number_stretch_level_start
822             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_start(i), ' m,'
823          ENDDO
824          WRITE( io, '(/,A)', advance='no') ' Vertical stretching starts at index: '
825          DO i = 1, number_stretch_level_start
826             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_start_index(i), ','
827          ENDDO
828          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at height:  '
829          DO i = 1, number_stretch_level_start
830             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_end(i), ' m,'
831          ENDDO
832          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at index:   '
833          DO i = 1, number_stretch_level_start
834             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_end_index(i), ','
835          ENDDO
836          WRITE( io, '(/,A)', advance='no') ' Factor used for stretching:          '
837          DO i = 1, number_stretch_level_start
838             WRITE ( io, '(F12.3,A1)', advance='no' )  dz_stretch_factor_array(i), ','
839          ENDDO
[94]840       ENDIF
[3065]841       
[94]842    ELSE
[3065]843       WRITE ( io, 250 )  dx, dy
844       DO i = 1, number_stretch_level_start+1
845          WRITE ( io, 253 )  i, dz(i)
846       ENDDO
847       
848       WRITE ( io, 251 ) (nx+1)*dx, (ny+1)*dy, zu(0)
849       
850       IF ( ANY( dz_stretch_level_start_index > 0 ) )  THEN
851          WRITE( io, '(A)', advance='no') ' Vertical stretching starts at height:'
852          DO i = 1, number_stretch_level_start
853             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_start(i), ' m,'
854          ENDDO
855          WRITE( io, '(/,A)', advance='no') ' Vertical stretching starts at index: '
856          DO i = 1, number_stretch_level_start
857             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_start_index(i), ','
858          ENDDO
859          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at height:  '
860          DO i = 1, number_stretch_level_start
861             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_end(i), ' m,'
862          ENDDO
863          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at index:   '
864          DO i = 1, number_stretch_level_start
865             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_end_index(i), ','
866          ENDDO
867          WRITE( io, '(/,A)', advance='no') ' Factor used for stretching:          '
868          DO i = 1, number_stretch_level_start
869             WRITE ( io, '(F12.3,A1)', advance='no' )  dz_stretch_factor_array(i), ','
870          ENDDO
[94]871       ENDIF
[1]872    ENDIF
[3065]873    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ),      &
[1]874                       MIN( nnz+2, nzt+2 )
875    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
876
877!
[1365]878!-- Large scale forcing and nudging
[2320]879    IF ( large_scale_forcing )  CALL lsf_nudging_header( io )
[1365]880
881!
882!-- Profile for the large scale vertial velocity
883!-- Building output strings, starting with surface value
884    IF ( large_scale_subsidence )  THEN
885       temperatures = '   0.0'
886       gradients = '------'
887       slices = '     0'
888       coordinates = '   0.0'
889       i = 1
890       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
891
892          WRITE (coor_chr,'(E10.2,7X)')  &
893                                w_subs(subs_vertical_gradient_level_i(i))
894          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
895
896          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
897          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
898
899          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
900          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
901
902          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
903          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
904
905          IF ( i == 10 )  THEN
906             EXIT
907          ELSE
908             i = i + 1
909          ENDIF
910
911       ENDDO
912
913 
914       IF ( .NOT. large_scale_forcing )  THEN
915          WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
916                             TRIM( gradients ), TRIM( slices )
917       ENDIF
918
919
920    ENDIF
921
922!-- Profile of the geostrophic wind (component ug)
923!-- Building output strings
924    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
925    gradients = '------'
926    slices = '     0'
927    coordinates = '   0.0'
928    i = 1
929    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
930     
931       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
932       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
933
934       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
935       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
936
937       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
938       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
939
940       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
941       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
942
943       IF ( i == 10 )  THEN
944          EXIT
945       ELSE
946          i = i + 1
947       ENDIF
948
949    ENDDO
950
951    IF ( .NOT. large_scale_forcing )  THEN
952       WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
953                          TRIM( gradients ), TRIM( slices )
954    ENDIF
955
956!-- Profile of the geostrophic wind (component vg)
957!-- Building output strings
958    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
959    gradients = '------'
960    slices = '     0'
961    coordinates = '   0.0'
962    i = 1
963    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
964
965       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
966       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
967
968       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
969       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
970
971       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
972       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
973
974       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
975       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
976
977       IF ( i == 10 )  THEN
978          EXIT
979       ELSE
980          i = i + 1
981       ENDIF
982 
983    ENDDO
984
985    IF ( .NOT. large_scale_forcing )  THEN
986       WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
987                          TRIM( gradients ), TRIM( slices )
988    ENDIF
989
990!
[1]991!-- Topography
992    WRITE ( io, 270 )  topography
993    SELECT CASE ( TRIM( topography ) )
994
995       CASE ( 'flat' )
996          ! no actions necessary
997
998       CASE ( 'single_building' )
999          blx = INT( building_length_x / dx )
1000          bly = INT( building_length_y / dy )
[1675]1001          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
1002          IF ( ABS( zw(bh  ) - building_height ) == &
1003               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
[1]1004
[1322]1005          IF ( building_wall_left == 9999999.9_wp )  THEN
[1]1006             building_wall_left = ( nx + 1 - blx ) / 2 * dx
1007          ENDIF
[1353]1008          bxl = INT ( building_wall_left / dx + 0.5_wp )
[1]1009          bxr = bxl + blx
1010
[1322]1011          IF ( building_wall_south == 9999999.9_wp )  THEN
[1]1012             building_wall_south = ( ny + 1 - bly ) / 2 * dy
1013          ENDIF
[1353]1014          bys = INT ( building_wall_south / dy + 0.5_wp )
[1]1015          byn = bys + bly
1016
1017          WRITE ( io, 271 )  building_length_x, building_length_y, &
1018                             building_height, bxl, bxr, bys, byn
1019
[240]1020       CASE ( 'single_street_canyon' )
[1675]1021          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
1022          IF ( ABS( zw(ch  ) - canyon_height ) == &
1023               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
[1322]1024          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[240]1025!
1026!--          Street canyon in y direction
1027             cwx = NINT( canyon_width_x / dx )
[1322]1028             IF ( canyon_wall_left == 9999999.9_wp )  THEN
[240]1029                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
1030             ENDIF
1031             cxl = NINT( canyon_wall_left / dx )
1032             cxr = cxl + cwx
1033             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
1034
[1322]1035          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[240]1036!
1037!--          Street canyon in x direction
1038             cwy = NINT( canyon_width_y / dy )
[1322]1039             IF ( canyon_wall_south == 9999999.9_wp )  THEN
[240]1040                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
1041             ENDIF
1042             cys = NINT( canyon_wall_south / dy )
1043             cyn = cys + cwy
1044             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
1045          ENDIF
1046
[2232]1047       CASE ( 'tunnel' )
1048          IF ( tunnel_width_x /= 9999999.9_wp )  THEN
1049!
1050!--          Tunnel axis in y direction
1051             IF ( tunnel_length == 9999999.9_wp  .OR.                          &
1052                  tunnel_length >= ( nx + 1 ) * dx )  THEN
1053                WRITE ( io, 273 )  'y', tunnel_height, tunnel_wall_depth,      &
1054                                        tunnel_width_x
1055             ELSE
1056                WRITE ( io, 274 )  'y', tunnel_height, tunnel_wall_depth,      &
1057                                        tunnel_width_x, tunnel_length
1058             ENDIF
1059
1060          ELSEIF ( tunnel_width_y /= 9999999.9_wp )  THEN
1061!
1062!--          Tunnel axis in x direction
1063             IF ( tunnel_length == 9999999.9_wp  .OR.                          &
1064                  tunnel_length >= ( ny + 1 ) * dy )  THEN
1065                WRITE ( io, 273 )  'x', tunnel_height, tunnel_wall_depth,      &
1066                                        tunnel_width_y
1067             ELSE
1068                WRITE ( io, 274 )  'x', tunnel_height, tunnel_wall_depth,      &
1069                                        tunnel_width_y, tunnel_length
1070             ENDIF
1071          ENDIF
1072
[1]1073    END SELECT
1074
[256]1075    IF ( TRIM( topography ) /= 'flat' )  THEN
1076       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
1077          IF ( TRIM( topography ) == 'single_building' .OR.  &
1078               TRIM( topography ) == 'single_street_canyon' )  THEN
1079             WRITE ( io, 278 )
1080          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
1081             WRITE ( io, 279 )
1082          ENDIF
1083       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
1084          WRITE ( io, 278 )
1085       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
1086          WRITE ( io, 279 )
1087       ENDIF
1088    ENDIF
1089
[2550]1090!-- Complex terrain
1091    IF ( complex_terrain )  THEN
1092       WRITE( io, 280 ) 
1093       IF ( turbulent_inflow )  THEN
[2698]1094          WRITE( io, 281 )  zu( get_topography_top_index_ji( 0, 0, 's' ) )
[2550]1095       ENDIF
1096       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1097          WRITE( io, 282 )
1098       ENDIF
1099    ENDIF
1100
[2776]1101    IF ( syn_turb_gen )  CALL stg_header ( io )
[2259]1102
[1826]1103    IF ( plant_canopy )  CALL pcm_header ( io )
[138]1104
[1817]1105    IF ( land_surface )  CALL lsm_header ( io )
[1484]1106
[1826]1107    IF ( radiation )  CALL radiation_header ( io )
[1551]1108
[2817]1109    IF ( gust_module_enabled )  CALL gust_header ( io )
1110
[1551]1111!
[1]1112!-- Boundary conditions
1113    IF ( ibc_p_b == 0 )  THEN
[1826]1114       r_lower = 'p(0)     = 0      |'
[1]1115    ELSEIF ( ibc_p_b == 1 )  THEN
[1826]1116       r_lower = 'p(0)     = p(1)   |'
[1]1117    ENDIF
1118    IF ( ibc_p_t == 0 )  THEN
[1826]1119       r_upper  = 'p(nzt+1) = 0      |'
[1]1120    ELSE
[1826]1121       r_upper  = 'p(nzt+1) = p(nzt) |'
[1]1122    ENDIF
1123
1124    IF ( ibc_uv_b == 0 )  THEN
[1826]1125       r_lower = TRIM( r_lower ) // ' uv(0)     = -uv(1)                |'
[1]1126    ELSE
[1826]1127       r_lower = TRIM( r_lower ) // ' uv(0)     = uv(1)                 |'
[1]1128    ENDIF
[132]1129    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
[1826]1130       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = 0                     |'
[132]1131    ELSEIF ( ibc_uv_t == 0 )  THEN
[1826]1132       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
[1]1133    ELSE
[1826]1134       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = uv(nzt)               |'
[1]1135    ENDIF
1136
1137    IF ( ibc_pt_b == 0 )  THEN
[1551]1138       IF ( land_surface )  THEN
[1826]1139          r_lower = TRIM( r_lower ) // ' pt(0)     = from soil model'
[1551]1140       ELSE
[1826]1141          r_lower = TRIM( r_lower ) // ' pt(0)     = pt_surface'
[1551]1142       ENDIF
[102]1143    ELSEIF ( ibc_pt_b == 1 )  THEN
[1826]1144       r_lower = TRIM( r_lower ) // ' pt(0)     = pt(1)'
[102]1145    ELSEIF ( ibc_pt_b == 2 )  THEN
[1826]1146       r_lower = TRIM( r_lower ) // ' pt(0)     = from coupled model'
[1]1147    ENDIF
1148    IF ( ibc_pt_t == 0 )  THEN
[1826]1149       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt_top'
[19]1150    ELSEIF( ibc_pt_t == 1 )  THEN
[1826]1151       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt)'
[19]1152    ELSEIF( ibc_pt_t == 2 )  THEN
[1826]1153       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
[667]1154
[1]1155    ENDIF
1156
[1826]1157    WRITE ( io, 300 )  r_lower, r_upper
[1]1158
1159    IF ( .NOT. constant_diffusion )  THEN
1160       IF ( ibc_e_b == 1 )  THEN
[1826]1161          r_lower = 'e(0)     = e(1)'
[1]1162       ELSE
[1826]1163          r_lower = 'e(0)     = e(1) = (u*/0.1)**2'
[1]1164       ENDIF
[1826]1165       r_upper = 'e(nzt+1) = e(nzt) = e(nzt-1)'
[1]1166
[1826]1167       WRITE ( io, 301 )  'e', r_lower, r_upper       
[1]1168
1169    ENDIF
1170
[97]1171    IF ( ocean )  THEN
[1826]1172       r_lower = 'sa(0)    = sa(1)'
[97]1173       IF ( ibc_sa_t == 0 )  THEN
[1826]1174          r_upper =  'sa(nzt+1) = sa_surface'
[1]1175       ELSE
[1826]1176          r_upper =  'sa(nzt+1) = sa(nzt)'
[1]1177       ENDIF
[1826]1178       WRITE ( io, 301 ) 'sa', r_lower, r_upper
[97]1179    ENDIF
[1]1180
[97]1181    IF ( humidity )  THEN
1182       IF ( ibc_q_b == 0 )  THEN
[1551]1183          IF ( land_surface )  THEN
[1826]1184             r_lower = 'q(0)     = from soil model'
[1551]1185          ELSE
[1826]1186             r_lower = 'q(0)     = q_surface'
[1551]1187          ENDIF
1188
[97]1189       ELSE
[1992]1190          r_lower = 'q(0)      = q(1)'
[97]1191       ENDIF
1192       IF ( ibc_q_t == 0 )  THEN
[1992]1193          r_upper =  'q(nzt+1) = q_top'
[97]1194       ELSE
[1992]1195          r_upper =  'q(nzt+1) = q(nzt) + dq/dz'
[97]1196       ENDIF
[1826]1197       WRITE ( io, 301 ) 'q', r_lower, r_upper
[97]1198    ENDIF
[1]1199
[97]1200    IF ( passive_scalar )  THEN
[1960]1201       IF ( ibc_s_b == 0 )  THEN
[1992]1202          r_lower = 's(0)      = s_surface'
[97]1203       ELSE
[1992]1204          r_lower = 's(0)      = s(1)'
[97]1205       ENDIF
[1960]1206       IF ( ibc_s_t == 0 )  THEN
[1992]1207          r_upper =  's(nzt+1) = s_top'
1208       ELSEIF ( ibc_s_t == 1 )  THEN
1209          r_upper =  's(nzt+1) = s(nzt)'
1210       ELSEIF ( ibc_s_t == 2 )  THEN
1211          r_upper =  's(nzt+1) = s(nzt) + ds/dz'
[97]1212       ENDIF
[1826]1213       WRITE ( io, 301 ) 's', r_lower, r_upper
[1]1214    ENDIF
1215
1216    IF ( use_surface_fluxes )  THEN
1217       WRITE ( io, 303 )
1218       IF ( constant_heatflux )  THEN
[1299]1219          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[2232]1220             IF ( surf_def_h(0)%ns >= 1 )  WRITE ( io, 306 )  surf_def_h(0)%shf(1)
[1241]1221          ELSE
1222             WRITE ( io, 306 )  surface_heatflux
1223          ENDIF
[1]1224          IF ( random_heatflux )  WRITE ( io, 307 )
1225       ENDIF
[75]1226       IF ( humidity  .AND.  constant_waterflux )  THEN
[1299]1227          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[2232]1228             WRITE ( io, 311 ) surf_def_h(0)%qsws(1)
[1241]1229          ELSE
1230             WRITE ( io, 311 ) surface_waterflux
1231          ENDIF
[1]1232       ENDIF
[1960]1233       IF ( passive_scalar  .AND.  constant_scalarflux )  THEN
1234          WRITE ( io, 313 ) surface_scalarflux
[1]1235       ENDIF
1236    ENDIF
1237
[19]1238    IF ( use_top_fluxes )  THEN
1239       WRITE ( io, 304 )
[102]1240       IF ( coupling_mode == 'uncoupled' )  THEN
[151]1241          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
[102]1242          IF ( constant_top_heatflux )  THEN
1243             WRITE ( io, 306 )  top_heatflux
1244          ENDIF
1245       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1246          WRITE ( io, 316 )
[19]1247       ENDIF
[1992]1248       IF ( ocean  .AND.  constant_top_salinityflux )                          &
[97]1249          WRITE ( io, 309 )  top_salinityflux
[1960]1250       IF ( humidity       )  WRITE ( io, 315 )
[1992]1251       IF ( passive_scalar .AND.  constant_top_scalarflux )                    &
1252          WRITE ( io, 302 ) top_scalarflux
[19]1253    ENDIF
1254
[1691]1255    IF ( constant_flux_layer )  THEN
1256       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length,                     &
1257                          z0h_factor*roughness_length, kappa,                  &
1258                          zeta_min, zeta_max
[1]1259       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
[75]1260       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
[1]1261          WRITE ( io, 312 )
1262       ENDIF
[1960]1263       IF ( passive_scalar  .AND.  .NOT. constant_scalarflux )  THEN
[1]1264          WRITE ( io, 314 )
1265       ENDIF
1266    ELSE
1267       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
[1691]1268          WRITE ( io, 310 )  zeta_min, zeta_max
[1]1269       ENDIF
1270    ENDIF
1271
1272    WRITE ( io, 317 )  bc_lr, bc_ns
[707]1273    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1159]1274       WRITE ( io, 318 )  use_cmax, pt_damping_width, pt_damping_factor       
[151]1275       IF ( turbulent_inflow )  THEN
[1560]1276          IF ( .NOT. recycling_yshift ) THEN
1277             WRITE ( io, 319 )  recycling_width, recycling_plane, &
1278                                inflow_damping_height, inflow_damping_width
1279          ELSE
1280             WRITE ( io, 322 )  recycling_width, recycling_plane, &
1281                                inflow_damping_height, inflow_damping_width
1282          END IF
[151]1283       ENDIF
[2050]1284       IF ( turbulent_outflow )  THEN
1285          WRITE ( io, 323 )  outflow_source_plane, INT(outflow_source_plane/dx)
1286       ENDIF
[1]1287    ENDIF
1288
1289!
[1365]1290!-- Initial Profiles
1291    WRITE ( io, 321 )
1292!
1293!-- Initial wind profiles
1294    IF ( u_profile(1) /= 9999999.9_wp )  WRITE ( io, 427 )
1295
1296!
1297!-- Initial temperature profile
1298!-- Building output strings, starting with surface temperature
1299    WRITE ( temperatures, '(F6.2)' )  pt_surface
1300    gradients = '------'
1301    slices = '     0'
1302    coordinates = '   0.0'
1303    i = 1
1304    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1305
1306       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1307       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1308
1309       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1310       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1311
1312       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1313       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1314
1315       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1316       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1317
1318       IF ( i == 10 )  THEN
1319          EXIT
1320       ELSE
1321          i = i + 1
1322       ENDIF
1323
1324    ENDDO
1325
1326    IF ( .NOT. nudging )  THEN
1327       WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1328                          TRIM( gradients ), TRIM( slices )
1329    ELSE
1330       WRITE ( io, 428 ) 
1331    ENDIF
1332
1333!
1334!-- Initial humidity profile
1335!-- Building output strings, starting with surface humidity
[1960]1336    IF ( humidity )  THEN
[1365]1337       WRITE ( temperatures, '(E8.1)' )  q_surface
1338       gradients = '--------'
1339       slices = '       0'
1340       coordinates = '     0.0'
1341       i = 1
1342       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1343         
1344          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1345          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1346
1347          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1348          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1349         
1350          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1351          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1352         
1353          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1354          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1355
1356          IF ( i == 10 )  THEN
1357             EXIT
1358          ELSE
1359             i = i + 1
1360          ENDIF
1361
1362       ENDDO
1363
[1960]1364       IF ( .NOT. nudging )  THEN
1365          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ),        &
[1365]1366                             TRIM( gradients ), TRIM( slices )
1367       ENDIF
1368    ENDIF
[1960]1369!
1370!-- Initial scalar profile
1371!-- Building output strings, starting with surface humidity
1372    IF ( passive_scalar )  THEN
1373       WRITE ( temperatures, '(E8.1)' )  s_surface
1374       gradients = '--------'
1375       slices = '       0'
1376       coordinates = '     0.0'
1377       i = 1
1378       DO  WHILE ( s_vertical_gradient_level_ind(i) /= -9999 )
1379         
[2073]1380          WRITE (coor_chr,'(E8.1,4X)')  s_init(s_vertical_gradient_level_ind(i))
[1960]1381          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
[1365]1382
[1960]1383          WRITE (coor_chr,'(E8.1,4X)')  s_vertical_gradient(i)
1384          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1385         
1386          WRITE (coor_chr,'(I8,4X)')  s_vertical_gradient_level_ind(i)
1387          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1388         
1389          WRITE (coor_chr,'(F8.1,4X)')  s_vertical_gradient_level(i)
1390          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1391
1392          IF ( i == 10 )  THEN
1393             EXIT
1394          ELSE
1395             i = i + 1
1396          ENDIF
1397
1398       ENDDO
1399
1400       WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ),           &
1401                          TRIM( gradients ), TRIM( slices )
1402    ENDIF   
1403
[1365]1404!
1405!-- Initial salinity profile
1406!-- Building output strings, starting with surface salinity
1407    IF ( ocean )  THEN
1408       WRITE ( temperatures, '(F6.2)' )  sa_surface
1409       gradients = '------'
1410       slices = '     0'
1411       coordinates = '   0.0'
1412       i = 1
1413       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1414
1415          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1416          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1417
1418          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1419          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1420
1421          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1422          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1423
1424          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1425          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1426
1427          IF ( i == 10 )  THEN
1428             EXIT
1429          ELSE
1430             i = i + 1
1431          ENDIF
1432
1433       ENDDO
1434
1435       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1436                          TRIM( gradients ), TRIM( slices )
1437    ENDIF
1438
1439
1440!
[1]1441!-- Listing of 1D-profiles
[151]1442    WRITE ( io, 325 )  dt_dopr_listing
[1353]1443    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1444       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1445    ENDIF
1446
1447!
1448!-- DATA output
1449    WRITE ( io, 330 )
[1353]1450    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1451       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1452    ENDIF
1453
1454!
1455!-- 1D-profiles
[346]1456    dopr_chr = 'Profile:'
[1]1457    IF ( dopr_n /= 0 )  THEN
1458       WRITE ( io, 331 )
1459
1460       output_format = ''
[1783]1461       output_format = netcdf_data_format_string
1462       IF ( netcdf_deflate == 0 )  THEN
1463          WRITE ( io, 344 )  output_format
1464       ELSE
1465          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1466       ENDIF
[1]1467
1468       DO  i = 1, dopr_n
1469          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
1470          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
1471             WRITE ( io, 332 )  dopr_chr
1472             dopr_chr = '       :'
1473          ENDIF
1474       ENDDO
1475
1476       IF ( dopr_chr /= '' )  THEN
1477          WRITE ( io, 332 )  dopr_chr
1478       ENDIF
1479       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
[1353]1480       IF ( skip_time_dopr /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dopr
[1]1481    ENDIF
1482
1483!
1484!-- 2D-arrays
1485    DO  av = 0, 1
1486
1487       i = 1
1488       do2d_xy = ''
1489       do2d_xz = ''
1490       do2d_yz = ''
1491       DO  WHILE ( do2d(av,i) /= ' ' )
1492
1493          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
1494          do2d_mode = do2d(av,i)(l-1:l)
1495
1496          SELECT CASE ( do2d_mode )
1497             CASE ( 'xy' )
1498                ll = LEN_TRIM( do2d_xy )
1499                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1500             CASE ( 'xz' )
1501                ll = LEN_TRIM( do2d_xz )
1502                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1503             CASE ( 'yz' )
1504                ll = LEN_TRIM( do2d_yz )
1505                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1506          END SELECT
1507
1508          i = i + 1
1509
1510       ENDDO
1511
1512       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
1513              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
[1327]1514              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) ) )  THEN
[1]1515
1516          IF (  av == 0 )  THEN
1517             WRITE ( io, 334 )  ''
1518          ELSE
1519             WRITE ( io, 334 )  '(time-averaged)'
1520          ENDIF
1521
1522          IF ( do2d_at_begin )  THEN
1523             begin_chr = 'and at the start'
1524          ELSE
1525             begin_chr = ''
1526          ENDIF
1527
1528          output_format = ''
[1783]1529          output_format = netcdf_data_format_string
1530          IF ( netcdf_deflate == 0 )  THEN
1531             WRITE ( io, 344 )  output_format
1532          ELSE
1533             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1534          ENDIF
[1]1535
1536          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
1537             i = 1
1538             slices = '/'
1539             coordinates = '/'
1540!
[1551]1541!--          Building strings with index and coordinate information of the
[1]1542!--          slices
1543             DO  WHILE ( section(i,1) /= -9999 )
1544
1545                WRITE (section_chr,'(I5)')  section(i,1)
1546                section_chr = ADJUSTL( section_chr )
1547                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1548
[206]1549                IF ( section(i,1) == -1 )  THEN
[1353]1550                   WRITE (coor_chr,'(F10.1)')  -1.0_wp
[206]1551                ELSE
1552                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
1553                ENDIF
[1]1554                coor_chr = ADJUSTL( coor_chr )
1555                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1556
1557                i = i + 1
1558             ENDDO
1559             IF ( av == 0 )  THEN
1560                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
1561                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
1562                                   TRIM( coordinates )
[1353]1563                IF ( skip_time_do2d_xy /= 0.0_wp )  THEN
[1]1564                   WRITE ( io, 339 )  skip_time_do2d_xy
1565                ENDIF
1566             ELSE
1567                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
1568                                   TRIM( begin_chr ), averaging_interval, &
1569                                   dt_averaging_input, 'k', TRIM( slices ), &
1570                                   TRIM( coordinates )
[1353]1571                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1572                   WRITE ( io, 339 )  skip_time_data_output_av
1573                ENDIF
1574             ENDIF
[1308]1575             IF ( netcdf_data_format > 4 )  THEN
1576                WRITE ( io, 352 )  ntdim_2d_xy(av)
1577             ELSE
1578                WRITE ( io, 353 )
1579             ENDIF
[1]1580          ENDIF
1581
1582          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
1583             i = 1
1584             slices = '/'
1585             coordinates = '/'
1586!
[1551]1587!--          Building strings with index and coordinate information of the
[1]1588!--          slices
1589             DO  WHILE ( section(i,2) /= -9999 )
1590
1591                WRITE (section_chr,'(I5)')  section(i,2)
1592                section_chr = ADJUSTL( section_chr )
1593                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1594
1595                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
1596                coor_chr = ADJUSTL( coor_chr )
1597                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1598
1599                i = i + 1
1600             ENDDO
1601             IF ( av == 0 )  THEN
1602                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
1603                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
1604                                   TRIM( coordinates )
[1353]1605                IF ( skip_time_do2d_xz /= 0.0_wp )  THEN
[1]1606                   WRITE ( io, 339 )  skip_time_do2d_xz
1607                ENDIF
1608             ELSE
1609                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
1610                                   TRIM( begin_chr ), averaging_interval, &
1611                                   dt_averaging_input, 'j', TRIM( slices ), &
1612                                   TRIM( coordinates )
[1353]1613                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1614                   WRITE ( io, 339 )  skip_time_data_output_av
1615                ENDIF
1616             ENDIF
[1308]1617             IF ( netcdf_data_format > 4 )  THEN
1618                WRITE ( io, 352 )  ntdim_2d_xz(av)
1619             ELSE
1620                WRITE ( io, 353 )
1621             ENDIF
[1]1622          ENDIF
1623
1624          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
1625             i = 1
1626             slices = '/'
1627             coordinates = '/'
1628!
[1551]1629!--          Building strings with index and coordinate information of the
[1]1630!--          slices
1631             DO  WHILE ( section(i,3) /= -9999 )
1632
1633                WRITE (section_chr,'(I5)')  section(i,3)
1634                section_chr = ADJUSTL( section_chr )
1635                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1636
1637                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
1638                coor_chr = ADJUSTL( coor_chr )
1639                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1640
1641                i = i + 1
1642             ENDDO
1643             IF ( av == 0 )  THEN
1644                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
1645                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
1646                                   TRIM( coordinates )
[1353]1647                IF ( skip_time_do2d_yz /= 0.0_wp )  THEN
[1]1648                   WRITE ( io, 339 )  skip_time_do2d_yz
1649                ENDIF
1650             ELSE
1651                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
1652                                   TRIM( begin_chr ), averaging_interval, &
1653                                   dt_averaging_input, 'i', TRIM( slices ), &
1654                                   TRIM( coordinates )
[1353]1655                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1656                   WRITE ( io, 339 )  skip_time_data_output_av
1657                ENDIF
1658             ENDIF
[1308]1659             IF ( netcdf_data_format > 4 )  THEN
1660                WRITE ( io, 352 )  ntdim_2d_yz(av)
1661             ELSE
1662                WRITE ( io, 353 )
1663             ENDIF
[1]1664          ENDIF
1665
1666       ENDIF
1667
1668    ENDDO
1669
1670!
1671!-- 3d-arrays
1672    DO  av = 0, 1
1673
1674       i = 1
1675       do3d_chr = ''
1676       DO  WHILE ( do3d(av,i) /= ' ' )
1677
1678          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
1679          i = i + 1
1680
1681       ENDDO
1682
1683       IF ( do3d_chr /= '' )  THEN
1684          IF ( av == 0 )  THEN
1685             WRITE ( io, 336 )  ''
1686          ELSE
1687             WRITE ( io, 336 )  '(time-averaged)'
1688          ENDIF
1689
[1783]1690          output_format = netcdf_data_format_string
1691          IF ( netcdf_deflate == 0 )  THEN
1692             WRITE ( io, 344 )  output_format
1693          ELSE
1694             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1695          ENDIF
[1]1696
1697          IF ( do3d_at_begin )  THEN
1698             begin_chr = 'and at the start'
1699          ELSE
1700             begin_chr = ''
1701          ENDIF
1702          IF ( av == 0 )  THEN
1703             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
1704                                zu(nz_do3d), nz_do3d
1705          ELSE
1706             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
1707                                TRIM( begin_chr ), averaging_interval, &
1708                                dt_averaging_input, zu(nz_do3d), nz_do3d
1709          ENDIF
1710
[1308]1711          IF ( netcdf_data_format > 4 )  THEN
1712             WRITE ( io, 352 )  ntdim_3d(av)
1713          ELSE
1714             WRITE ( io, 353 )
1715          ENDIF
1716
[1]1717          IF ( av == 0 )  THEN
[1353]1718             IF ( skip_time_do3d /= 0.0_wp )  THEN
[1]1719                WRITE ( io, 339 )  skip_time_do3d
1720             ENDIF
1721          ELSE
[1353]1722             IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1723                WRITE ( io, 339 )  skip_time_data_output_av
1724             ENDIF
1725          ENDIF
1726
1727       ENDIF
1728
1729    ENDDO
1730
1731!
[410]1732!-- masked arrays
1733    IF ( masks > 0 )  WRITE ( io, 345 )  &
1734         mask_scale_x, mask_scale_y, mask_scale_z
1735    DO  mid = 1, masks
1736       DO  av = 0, 1
1737
1738          i = 1
1739          domask_chr = ''
1740          DO  WHILE ( domask(mid,av,i) /= ' ' )
1741             domask_chr = TRIM( domask_chr ) // ' ' //  &
1742                          TRIM( domask(mid,av,i) ) // ','
1743             i = i + 1
1744          ENDDO
1745
1746          IF ( domask_chr /= '' )  THEN
1747             IF ( av == 0 )  THEN
1748                WRITE ( io, 346 )  '', mid
1749             ELSE
1750                WRITE ( io, 346 )  ' (time-averaged)', mid
1751             ENDIF
1752
[1783]1753             output_format = netcdf_data_format_string
[1308]1754!--          Parallel output not implemented for mask data, hence
1755!--          output_format must be adjusted.
1756             IF ( netcdf_data_format == 5 ) output_format = 'netCDF4/HDF5'
1757             IF ( netcdf_data_format == 6 ) output_format = 'netCDF4/HDF5 classic'
[1783]1758             IF ( netcdf_deflate == 0 )  THEN
1759                WRITE ( io, 344 )  output_format
1760             ELSE
1761                WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1762             ENDIF
[410]1763
1764             IF ( av == 0 )  THEN
1765                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1766             ELSE
1767                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1768                                   averaging_interval, dt_averaging_input
1769             ENDIF
1770
1771             IF ( av == 0 )  THEN
[1353]1772                IF ( skip_time_domask(mid) /= 0.0_wp )  THEN
[410]1773                   WRITE ( io, 339 )  skip_time_domask(mid)
1774                ENDIF
1775             ELSE
[1353]1776                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[410]1777                   WRITE ( io, 339 )  skip_time_data_output_av
1778                ENDIF
1779             ENDIF
1780!
1781!--          output locations
1782             DO  dim = 1, 3
[1353]1783                IF ( mask(mid,dim,1) >= 0.0_wp )  THEN
[410]1784                   count = 0
[1353]1785                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0_wp )
[410]1786                      count = count + 1
1787                   ENDDO
1788                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1789                                      mask(mid,dim,:count)
[1353]1790                ELSEIF ( mask_loop(mid,dim,1) < 0.0_wp .AND.  &
1791                         mask_loop(mid,dim,2) < 0.0_wp .AND.  &
1792                         mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1793                   WRITE ( io, 350 )  dir(dim), dir(dim)
[1353]1794                ELSEIF ( mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1795                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1796                                      mask_loop(mid,dim,1:2)
1797                ELSE
1798                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1799                                      mask_loop(mid,dim,1:3)
1800                ENDIF
1801             ENDDO
1802          ENDIF
1803
1804       ENDDO
1805    ENDDO
1806
1807!
[1]1808!-- Timeseries
[1322]1809    IF ( dt_dots /= 9999999.9_wp )  THEN
[1]1810       WRITE ( io, 340 )
1811
[1783]1812       output_format = netcdf_data_format_string
1813       IF ( netcdf_deflate == 0 )  THEN
1814          WRITE ( io, 344 )  output_format
1815       ELSE
1816          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1817       ENDIF
[1]1818       WRITE ( io, 341 )  dt_dots
1819    ENDIF
1820
1821#if defined( __dvrp_graphics )
1822!
1823!-- Dvrp-output
[1322]1824    IF ( dt_dvrp /= 9999999.9_wp )  THEN
[1]1825       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1826                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1827       i = 1
1828       l = 0
[336]1829       m = 0
[1]1830       DO WHILE ( mode_dvrp(i) /= ' ' )
1831          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
[130]1832             READ ( mode_dvrp(i), '(10X,I2)' )  j
[1]1833             l = l + 1
1834             IF ( do3d(0,j) /= ' ' )  THEN
[336]1835                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1836                                   isosurface_color(:,l)
[1]1837             ENDIF
1838          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
[130]1839             READ ( mode_dvrp(i), '(6X,I2)' )  j
[336]1840             m = m + 1
1841             IF ( do2d(0,j) /= ' ' )  THEN
1842                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1843                                   slicer_range_limits_dvrp(:,m)
1844             ENDIF
[1]1845          ENDIF
1846          i = i + 1
1847       ENDDO
[237]1848
[336]1849       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1850                          superelevation_y, superelevation, clip_dvrp_l, &
1851                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1852
1853       IF ( TRIM( topography ) /= 'flat' )  THEN
1854          WRITE ( io, 366 )  topography_color
1855          IF ( cluster_size > 1 )  THEN
1856             WRITE ( io, 367 )  cluster_size
1857          ENDIF
[237]1858       ENDIF
1859
[1]1860    ENDIF
1861#endif
[1957]1862!
1863!-- Output of virtual flight information
1864    IF ( virtual_flight )  CALL flight_header( io )
[1]1865
1866!
[1833]1867!-- Output of spectra related quantities
1868    IF ( calculate_spectra )  CALL spectra_header( io )
[1]1869
1870    WRITE ( io, 99 )
1871
1872!
1873!-- Physical quantities
1874    WRITE ( io, 400 )
1875
1876!
1877!-- Geostrophic parameters
[2575]1878    WRITE ( io, 410 )  latitude, longitude, omega, f, fs
[1]1879
[2544]1880 !
1881!-- Geostrophic parameters
1882    WRITE ( io, 456 )  day_of_year_init, time_utc_init
1883   
[1]1884!
1885!-- Other quantities
1886    WRITE ( io, 411 )  g
[1551]1887
[1179]1888    WRITE ( io, 412 )  TRIM( reference_state )
1889    IF ( use_single_reference_value )  THEN
[97]1890       IF ( ocean )  THEN
[1179]1891          WRITE ( io, 413 )  prho_reference
[97]1892       ELSE
[1179]1893          WRITE ( io, 414 )  pt_reference
[97]1894       ENDIF
1895    ENDIF
[1]1896
1897!
1898!-- Cloud physics parameters
[1299]1899    IF ( cloud_physics )  THEN
[57]1900       WRITE ( io, 415 )
1901       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
[1822]1902       IF ( microphysics_seifert )  THEN
[1353]1903          WRITE ( io, 510 ) 1.0E-6_wp * nc_const
[1822]1904          WRITE ( io, 511 ) c_sedimentation
[1115]1905       ENDIF
[1]1906    ENDIF
1907
1908!
[824]1909!-- Cloud physcis parameters / quantities / numerical methods
1910    WRITE ( io, 430 )
1911    IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets)  THEN
1912       WRITE ( io, 431 )
1913    ELSEIF ( humidity  .AND.  cloud_physics )  THEN
1914       WRITE ( io, 432 )
[1496]1915       IF ( cloud_top_radiation )  WRITE ( io, 132 )
[1822]1916       IF ( microphysics_kessler )  THEN
1917          WRITE ( io, 133 )
1918       ELSEIF ( microphysics_seifert )  THEN
[1831]1919          IF ( cloud_water_sedimentation )  WRITE ( io, 506 )
[1822]1920          WRITE ( io, 505 )
[1831]1921          IF ( collision_turbulence )  WRITE ( io, 507 )
[1822]1922          IF ( ventilation_effect )  WRITE ( io, 508 )
1923          IF ( limiter_sedimentation )  WRITE ( io, 509 )
[1115]1924       ENDIF
[824]1925    ELSEIF ( humidity  .AND.  cloud_droplets )  THEN
1926       WRITE ( io, 433 )
1927       IF ( curvature_solution_effects )  WRITE ( io, 434 )
[825]1928       IF ( collision_kernel /= 'none' )  THEN
1929          WRITE ( io, 435 )  TRIM( collision_kernel )
[828]1930          IF ( collision_kernel(6:9) == 'fast' )  THEN
1931             WRITE ( io, 436 )  radius_classes, dissipation_classes
1932          ENDIF
[825]1933       ELSE
[828]1934          WRITE ( io, 437 )
[825]1935       ENDIF
[824]1936    ENDIF
1937
1938!
[1]1939!-- LES / turbulence parameters
1940    WRITE ( io, 450 )
1941
1942!--
1943! ... LES-constants used must still be added here
1944!--
1945    IF ( constant_diffusion )  THEN
1946       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1947                          prandtl_number
1948    ENDIF
1949    IF ( .NOT. constant_diffusion)  THEN
[1353]1950       IF ( e_init > 0.0_wp )  WRITE ( io, 455 )  e_init
1951       IF ( e_min > 0.0_wp )  WRITE ( io, 454 )  e_min
[1]1952       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1953    ENDIF
[3083]1954    IF ( rans_mode )  THEN
1955       WRITE ( io, 457 )  rans_const_c, rans_const_sigma
1956    ENDIF
[1]1957!
1958!-- Special actions during the run
1959    WRITE ( io, 470 )
1960    IF ( create_disturbances )  THEN
1961       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1962                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1963                          zu(disturbance_level_ind_t), disturbance_level_ind_t
[707]1964       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1]1965          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1966       ELSE
1967          WRITE ( io, 473 )  disturbance_energy_limit
1968       ENDIF
1969       WRITE ( io, 474 )  TRIM( random_generator )
1970    ENDIF
[1353]1971    IF ( pt_surface_initial_change /= 0.0_wp )  THEN
[1]1972       WRITE ( io, 475 )  pt_surface_initial_change
1973    ENDIF
[1353]1974    IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1975       WRITE ( io, 476 )  q_surface_initial_change       
1976    ENDIF
[1353]1977    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1978       WRITE ( io, 477 )  q_surface_initial_change       
1979    ENDIF
1980
[60]1981    IF ( particle_advection )  THEN
[1]1982!
[60]1983!--    Particle attributes
1984       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1985                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
[1359]1986                          end_time_prel
[60]1987       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1988       IF ( random_start_position )  WRITE ( io, 481 )
[1575]1989       IF ( seed_follows_topography )  WRITE ( io, 496 )
[60]1990       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1991       WRITE ( io, 495 )  total_number_of_particles
[1322]1992       IF ( dt_write_particle_data /= 9999999.9_wp )  THEN
[60]1993          WRITE ( io, 485 )  dt_write_particle_data
[1327]1994          IF ( netcdf_data_format > 1 )  THEN
1995             output_format = 'netcdf (64 bit offset) and binary'
[1]1996          ELSE
[1327]1997             output_format = 'netcdf and binary'
[1]1998          ENDIF
[1783]1999          IF ( netcdf_deflate == 0 )  THEN
2000             WRITE ( io, 344 )  output_format
2001          ELSE
2002             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
2003          ENDIF
[1]2004       ENDIF
[1322]2005       IF ( dt_dopts /= 9999999.9_wp )  WRITE ( io, 494 )  dt_dopts
[60]2006       IF ( write_particle_statistics )  WRITE ( io, 486 )
[1]2007
[60]2008       WRITE ( io, 487 )  number_of_particle_groups
[1]2009
[60]2010       DO  i = 1, number_of_particle_groups
[1322]2011          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9_wp )  THEN
[1353]2012             WRITE ( io, 490 )  i, 0.0_wp
[60]2013             WRITE ( io, 492 )
[1]2014          ELSE
[60]2015             WRITE ( io, 490 )  i, radius(i)
[1353]2016             IF ( density_ratio(i) /= 0.0_wp )  THEN
[60]2017                WRITE ( io, 491 )  density_ratio(i)
2018             ELSE
2019                WRITE ( io, 492 )
2020             ENDIF
[1]2021          ENDIF
[60]2022          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
2023                             pdx(i), pdy(i), pdz(i)
[336]2024          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
[60]2025       ENDDO
[1]2026
[60]2027    ENDIF
[1]2028
[60]2029
[1]2030!
2031!-- Parameters of 1D-model
2032    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
2033       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
2034                          mixing_length_1d, dissipation_1d
2035       IF ( damp_level_ind_1d /= nzt+1 )  THEN
2036          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
2037       ENDIF
2038    ENDIF
2039
2040!
[1551]2041!-- User-defined information
[1]2042    CALL user_header( io )
2043
2044    WRITE ( io, 99 )
2045
2046!
2047!-- Write buffer contents to disc immediately
[1808]2048    FLUSH( io )
[1]2049
2050!
2051!-- Here the FORMATs start
2052
2053 99 FORMAT (1X,78('-'))
[1468]2054100 FORMAT (/1X,'******************************',4X,44('-')/        &
2055            1X,'* ',A,' *',4X,A/                               &
2056            1X,'******************************',4X,44('-'))
[2298]2057101 FORMAT (35X,'coupled run: ',A/ &
[1468]2058            35X,42('-'))
2059102 FORMAT (/' Date:                 ',A8,4X,'Run:       ',A20/      &
2060            ' Time:                 ',A8,4X,'Run-No.:   ',I2.2/     &
[1106]2061            ' Run on host:        ',A10)
[1]2062#if defined( __parallel )
[1468]2063103 FORMAT (' Number of PEs:',10X,I6,4X,'Processor grid (x,y): (',I4,',',I4, &
[1]2064              ')',1X,A)
[1468]2065104 FORMAT (' Number of PEs:',10X,I6,4X,'Tasks:',I4,'   threads per task:',I4/ &
2066              35X,'Processor grid (x,y): (',I4,',',I4,')',1X,A)
2067105 FORMAT (35X,'One additional PE is used to handle'/37X,'the dvrp output!')
2068107 FORMAT (35X,'A 1d-decomposition along ',A,' is used')
2069108 FORMAT (35X,'Max. # of parallel I/O streams is ',I5)
2070109 FORMAT (35X,'Precursor run for coupled atmos-ocean run'/ &
2071            35X,42('-'))
2072114 FORMAT (35X,'Coupled atmosphere-ocean run following'/ &
2073            35X,'independent precursor runs'/             &
2074            35X,42('-'))
[1]2075#endif
2076110 FORMAT (/' Numerical Schemes:'/ &
2077             ' -----------------'/)
[2696]2078124 FORMAT (' --> Use the ',A,' turbulence closure (',A,' mode).')
[2037]2079121 FORMAT (' --> Use the ',A,' approximation for the model equations.')
[1]2080111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
2081112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
[1697]2082            '     Iterations (initial/other): ',I3,'/',I3,'  omega =',F6.3)
[1]2083113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
2084                  ' or Upstream')
[1216]2085115 FORMAT ('     FFT and transpositions are overlapping')
[1]2086116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
2087                  ' or Upstream')
2088118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
[1106]2089119 FORMAT (' --> Galilei-Transform applied to horizontal advection:'/ &
2090            '     translation velocity = ',A/ &
[1]2091            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
2092122 FORMAT (' --> Time differencing scheme: ',A)
[108]2093123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
[1697]2094            '     maximum damping coefficient:',F6.3, ' 1/s')
[1]2095129 FORMAT (' --> Additional prognostic equation for the specific humidity')
2096130 FORMAT (' --> Additional prognostic equation for the total water content')
[940]2097131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
2098                  F6.2, ' K assumed')
[824]2099132 FORMAT ('     Parameterization of long-wave radiation processes via'/ &
[1]2100            '     effective emissivity scheme')
[824]2101133 FORMAT ('     Precipitation parameterization via Kessler-Scheme')
[1]2102134 FORMAT (' --> Additional prognostic equation for a passive scalar')
[1575]2103135 FORMAT (' --> Solve perturbation pressure via ',A,' method (', &
[1]2104                  A,'-cycle)'/ &
2105            '     number of grid levels:                   ',I2/ &
2106            '     Gauss-Seidel red/black iterations:       ',I2)
2107136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
2108                  I3,')')
2109137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
2110            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
2111                  I3,')'/ &
2112            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
2113                  I3,')')
[63]2114139 FORMAT (' --> Loop optimization method: ',A)
[1]2115140 FORMAT ('     maximum residual allowed:                ',E10.3)
2116141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
2117142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
2118                  'step')
[87]2119143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
2120                  'kinetic energy')
[927]2121144 FORMAT ('     masking method is used')
[1]2122150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
[241]2123                  'conserved'/ &
2124            '     using the ',A,' mode')
2125151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
[306]2126152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
2127           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
2128           /'     starting from dp_level_b =', F8.3, 'm', A /)
[1]2129200 FORMAT (//' Run time and time step information:'/ &
2130             ' ----------------------------------'/)
[1106]2131201 FORMAT ( ' Timestep:             variable     maximum value: ',F6.3,' s', &
[1697]2132             '    CFL-factor:',F5.2)
[1106]2133202 FORMAT ( ' Timestep:          dt = ',F6.3,' s'/)
2134203 FORMAT ( ' Start time:          ',F9.3,' s'/ &
2135             ' End time:            ',F9.3,' s')
[1]2136204 FORMAT ( A,F9.3,' s')
2137205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
[1106]2138206 FORMAT (/' Time reached:        ',F9.3,' s'/ &
2139             ' CPU-time used:       ',F9.3,' s     per timestep:               ', &
2140               '  ',F9.3,' s'/                                                    &
[1111]2141             '                                      per second of simulated tim', &
[1]2142               'e: ',F9.3,' s')
[1106]2143207 FORMAT ( ' Coupling start time: ',F9.3,' s')
[1]2144250 FORMAT (//' Computational grid and domain size:'/ &
2145              ' ----------------------------------'// &
[3065]2146              ' Grid length:      dx =    ',F8.3,' m    dy =    ',F8.3, ' m')
2147251 FORMAT (  /' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
[1]2148              ' m  z(u) = ',F10.3,' m'/)
[3065]2149253 FORMAT ( '                dz(',I1,') =    ', F8.3, ' m')
2150254 FORMAT (//' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
[1]2151            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
2152260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
2153             ' degrees')
[1551]2154270 FORMAT (//' Topography information:'/ &
2155              ' ----------------------'// &
[1]2156              1X,'Topography: ',A)
2157271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
2158              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
2159                ' / ',I4)
[240]2160272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
2161              ' direction' / &
2162              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
2163              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
[2232]2164273 FORMAT (  ' Tunnel of infinite length in ',A, &
2165              ' direction' / &
2166              ' Tunnel height: ', F6.2, / &
2167              ' Tunnel-wall depth: ', F6.2      / &
2168              ' Tunnel width: ', F6.2 )
2169274 FORMAT (  ' Tunnel in ', A, ' direction.' / &
2170              ' Tunnel height: ', F6.2, / &   
2171              ' Tunnel-wall depth: ', F6.2      / &
2172              ' Tunnel width: ', F6.2, / &
2173              ' Tunnel length: ', F6.2 )
[256]2174278 FORMAT (' Topography grid definition convention:'/ &
2175            ' cell edge (staggered grid points'/  &
2176            ' (u in x-direction, v in y-direction))' /)
2177279 FORMAT (' Topography grid definition convention:'/ &
2178            ' cell center (scalar grid points)' /)
[2550]2179280 FORMAT (' Complex terrain simulation is activated.')
2180281 FORMAT ('    --> Mean inflow profiles are adjusted.' / &
2181            '    --> Elevation of inflow boundary: ', F7.1, ' m' )
2182282 FORMAT ('    --> Initial data from 3D-precursor run is shifted' / &
2183            '        vertically depending on local surface height.')
[1]2184300 FORMAT (//' Boundary conditions:'/ &
2185             ' -------------------'// &
2186             '                     p                    uv             ', &
[1551]2187             '                     pt'// &
[1]2188             ' B. bound.: ',A/ &
2189             ' T. bound.: ',A)
[97]2190301 FORMAT (/'                     ',A// &
[1]2191             ' B. bound.: ',A/ &
2192             ' T. bound.: ',A)
[19]2193303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
2194304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
[2270]2195305 FORMAT (//'    Constant flux layer between bottom surface and first ',     &
2196              'computational u,v-level:'// &
2197             '       z_mo = ',F6.2,' m   z0 =',F7.4,' m   z0h =',F8.5,&
[1697]2198             ' m   kappa =',F5.2/ &
2199             '       Rif value range:   ',F8.2,' <= rif <=',F6.2)
[97]2200306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
[1]2201307 FORMAT ('       Heatflux has a random normal distribution')
2202308 FORMAT ('       Predefined surface temperature')
[97]2203309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
[1]2204310 FORMAT (//'    1D-Model:'// &
2205             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
[1960]2206311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' kg/kg m/s')
[1]2207312 FORMAT ('       Predefined surface humidity')
2208313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
2209314 FORMAT ('       Predefined scalar value at the surface')
[1992]2210302 FORMAT ('       Predefined constant scalarflux:   ',F9.6,' kg/(m**2 s)')
2211315 FORMAT ('       Humidity flux at top surface is 0.0')
[102]2212316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
2213                    'atmosphere model')
[1]2214317 FORMAT (//' Lateral boundaries:'/ &
2215            '       left/right:  ',A/    &
2216            '       north/south: ',A)
[1159]2217318 FORMAT (/'       use_cmax: ',L1 / &
2218            '       pt damping layer width = ',F8.2,' m, pt ', &
[1697]2219                    'damping factor =',F7.4)
[151]2220319 FORMAT ('       turbulence recycling at inflow switched on'/ &
2221            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
2222            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
2223320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
[103]2224            '                                          v: ',F9.6,' m**2/s**2')
[1365]2225321 FORMAT (//' Initial profiles:'/ &
2226              ' ----------------')
[1560]2227322 FORMAT ('       turbulence recycling at inflow switched on'/ &
2228            '       y shift of the recycled inflow turbulence switched on'/ &
2229            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
[1592]2230            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m'/)
[2050]2231323 FORMAT ('       turbulent outflow conditon switched on'/ &
2232            '       position of outflow source plane: ',F7.1,' m   ', &
2233                    'grid index: ', I4)
[151]2234325 FORMAT (//' List output:'/ &
[1]2235             ' -----------'//  &
2236            '    1D-Profiles:'/    &
[2883]2237            '       Output every             ',F10.2,' s')
[151]2238326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
[1]2239            '       Averaging input every    ',F8.2,' s')
2240330 FORMAT (//' Data output:'/ &
2241             ' -----------'/)
2242331 FORMAT (/'    1D-Profiles:')
2243332 FORMAT (/'       ',A)
2244333 FORMAT ('       Output every             ',F8.2,' s',/ &
2245            '       Time averaged over       ',F8.2,' s'/ &
2246            '       Averaging input every    ',F8.2,' s')
2247334 FORMAT (/'    2D-Arrays',A,':')
2248335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2249            '       Output every             ',F8.2,' s  ',A/ &
2250            '       Cross sections at ',A1,' = ',A/ &
2251            '       scalar-coordinates:   ',A,' m'/)
2252336 FORMAT (/'    3D-Arrays',A,':')
2253337 FORMAT (/'       Arrays: ',A/ &
2254            '       Output every             ',F8.2,' s  ',A/ &
2255            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
2256339 FORMAT ('       No output during initial ',F8.2,' s')
2257340 FORMAT (/'    Time series:')
2258341 FORMAT ('       Output every             ',F8.2,' s'/)
2259342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2260            '       Output every             ',F8.2,' s  ',A/ &
2261            '       Time averaged over       ',F8.2,' s'/ &
2262            '       Averaging input every    ',F8.2,' s'/ &
2263            '       Cross sections at ',A1,' = ',A/ &
2264            '       scalar-coordinates:   ',A,' m'/)
2265343 FORMAT (/'       Arrays: ',A/ &
2266            '       Output every             ',F8.2,' s  ',A/ &
2267            '       Time averaged over       ',F8.2,' s'/ &
2268            '       Averaging input every    ',F8.2,' s'/ &
2269            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
[292]2270344 FORMAT ('       Output format: ',A/)
[410]2271345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
2272            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
2273            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
2274            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
2275346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
2276347 FORMAT ('       Variables: ',A/ &
2277            '       Output every             ',F8.2,' s')
2278348 FORMAT ('       Variables: ',A/ &
2279            '       Output every             ',F8.2,' s'/ &
2280            '       Time averaged over       ',F8.2,' s'/ &
2281            '       Averaging input every    ',F8.2,' s')
2282349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2283            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
2284            13('       ',8(F8.2,',')/) )
2285350 FORMAT (/'       Output locations in ',A,'-direction: ', &
2286            'all gridpoints along ',A,'-direction (default).' )
2287351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2288            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
2289            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
[1313]2290352 FORMAT  (/'       Number of output time levels allowed: ',I3 /)
2291353 FORMAT  (/'       Number of output time levels allowed: unlimited' /)
[1783]2292354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
[1]2293#if defined( __dvrp_graphics )
2294360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
2295            '       Output every      ',F7.1,' s'/ &
2296            '       Output mode:      ',A/ &
2297            '       Host / User:      ',A,' / ',A/ &
2298            '       Directory:        ',A// &
2299            '       The sequence contains:')
[337]2300361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
2301            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
2302362 FORMAT (/'       Slicer plane ',A/ &
[336]2303            '       Slicer limits: [',F6.2,',',F6.2,']')
[337]2304365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
[336]2305            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
2306                     ')'/ &
2307            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
2308            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
[337]2309366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
[336]2310367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
[1]2311#endif
2312400 FORMAT (//' Physical quantities:'/ &
2313              ' -------------------'/)
[2575]2314410 FORMAT ('    Geograph. latitude  :   latitude  = ',F4.1,' degr'/   &
2315            '    Geograph. longitude :   longitude = ',F4.1,' degr'/   &
[1697]2316            '    Angular velocity    :   omega  =',E10.3,' rad/s'/  &
[1551]2317            '    Coriolis parameter  :   f      = ',F9.6,' 1/s'/    &
2318            '                            f*     = ',F9.6,' 1/s')
2319411 FORMAT (/'    Gravity             :   g      = ',F4.1,' m/s**2')
[1179]2320412 FORMAT (/'    Reference state used in buoyancy terms: ',A)
2321413 FORMAT ('       Reference density in buoyancy terms: ',F8.3,' kg/m**3')
2322414 FORMAT ('       Reference temperature in buoyancy terms: ',F8.4,' K')
[1551]2323415 FORMAT (/' Cloud physics parameters:'/ &
2324             ' ------------------------'/)
2325416 FORMAT ('    Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
2326            '    Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
[1697]2327            '    Density of air     :   rho_0 =',F6.3,' kg/m**3'/  &
[1551]2328            '    Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
[1697]2329            '    Vapourization heat :   L_v   =',E9.2,' J/kg')
[1551]2330418 FORMAT (/'    Day of the year at model start :   day_init      =     ',I3 &
2331            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
[1]2332420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
2333            '       Height:        ',A,'  m'/ &
2334            '       Temperature:   ',A,'  K'/ &
2335            '       Gradient:      ',A,'  K/100m'/ &
2336            '       Gridpoint:     ',A)
2337421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
2338            '       Height:      ',A,'  m'/ &
2339            '       Humidity:    ',A,'  kg/kg'/ &
2340            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
2341            '       Gridpoint:   ',A)
2342422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
2343            '       Height:                  ',A,'  m'/ &
2344            '       Scalar concentration:    ',A,'  kg/m**3'/ &
2345            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
2346            '       Gridpoint:               ',A)
2347423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
2348            '       Height:      ',A,'  m'/ &
2349            '       ug:          ',A,'  m/s'/ &
2350            '       Gradient:    ',A,'  1/100s'/ &
2351            '       Gridpoint:   ',A)
2352424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
2353            '       Height:      ',A,'  m'/ &
[97]2354            '       vg:          ',A,'  m/s'/ &
[1]2355            '       Gradient:    ',A,'  1/100s'/ &
2356            '       Gridpoint:   ',A)
[97]2357425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
2358            '       Height:     ',A,'  m'/ &
2359            '       Salinity:   ',A,'  psu'/ &
2360            '       Gradient:   ',A,'  psu/100m'/ &
2361            '       Gridpoint:  ',A)
[411]2362426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
2363            '       Height:      ',A,'  m'/ &
2364            '       w_subs:      ',A,'  m/s'/ &
2365            '       Gradient:    ',A,'  (m/s)/100m'/ &
2366            '       Gridpoint:   ',A)
[767]2367427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
2368                  ' profiles')
[1241]2369428 FORMAT (/'    Initial profiles (u, v, pt, q) are taken from file '/ &
2370             '    NUDGING_DATA')
[824]2371430 FORMAT (//' Cloud physics quantities / methods:'/ &
2372              ' ----------------------------------'/)
[1960]2373431 FORMAT ('    Humidity is considered, bu no condensation')
[824]2374432 FORMAT ('    Bulk scheme with liquid water potential temperature and'/ &
2375            '    total water content is used.'/ &
2376            '    Condensation is parameterized via 0% - or 100% scheme.')
2377433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
2378                 'icle model')
2379434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
2380                 ' droplets < 1.0E-6 m')
[825]2381435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
[828]2382436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
2383                    'are used'/ &
2384            '          number of radius classes:       ',I3,'    interval ', &
2385                       '[1.0E-6,2.0E-4] m'/ &
2386            '          number of dissipation classes:   ',I2,'    interval ', &
2387                       '[0,1000] cm**2/s**3')
2388437 FORMAT ('    Droplet collision is switched off')
[1]2389450 FORMAT (//' LES / Turbulence quantities:'/ &
2390              ' ---------------------------'/)
[824]2391451 FORMAT ('    Diffusion coefficients are constant:'/ &
2392            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
[1697]2393453 FORMAT ('    Mixing length is limited to',F5.2,' * z')
[824]2394454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
2395455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
[2544]2396456 FORMAT ('    Day of the year at model start :   day_init = ',I3             &
2397            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
[3083]2398457 FORMAT ('    RANS-mode constants: c_0 = ',F9.5/         &
2399            '                         c_1 = ',F9.5/         &
2400            '                         c_2 = ',F9.5/         &
2401            '                         c_3 = ',F9.5/         &
2402            '                         c_4 = ',F9.5/         &
2403            '                         sigma_e    = ',F9.5/  &
2404            '                         sigma_diss = ',F9.5)
[1]2405470 FORMAT (//' Actions during the simulation:'/ &
2406              ' -----------------------------'/)
[94]2407471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
[1697]2408            '    Disturbance amplitude           :    ',F5.2, ' m/s'/       &
[94]2409            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
2410            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
[1]2411472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
2412                 ' to i/j =',I4)
2413473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
[1697]2414                 F6.3, ' m**2/s**2')
[1]2415474 FORMAT ('    Random number generator used    : ',A/)
2416475 FORMAT ('    The surface temperature is increased (or decreased, ', &
2417                 'respectively, if'/ &
2418            '    the value is negative) by ',F5.2,' K at the beginning of the',&
2419                 ' 3D-simulation'/)
2420476 FORMAT ('    The surface humidity is increased (or decreased, ',&
2421                 'respectively, if the'/ &
2422            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
2423                 ' the 3D-simulation'/)
2424477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
2425                 'respectively, if the'/ &
2426            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
2427                 ' the 3D-simulation'/)
2428480 FORMAT ('    Particles:'/ &
2429            '    ---------'// &
2430            '       Particle advection is active (switched on at t = ', F7.1, &
2431                    ' s)'/ &
2432            '       Start of new particle generations every  ',F6.1,' s'/ &
2433            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
2434            '                            bottom:     ', A, ' top:         ', A/&
2435            '       Maximum particle age:                 ',F9.1,' s'/ &
[1359]2436            '       Advection stopped at t = ',F9.1,' s'/)
[1]2437481 FORMAT ('       Particles have random start positions'/)
[336]2438482 FORMAT ('          Particles are advected only horizontally'/)
[1]2439485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
2440486 FORMAT ('       Particle statistics are written on file'/)
2441487 FORMAT ('       Number of particle groups: ',I2/)
2442488 FORMAT ('       SGS velocity components are used for particle advection'/ &
[1697]2443            '          minimum timestep for advection:', F8.5/)
[1]2444489 FORMAT ('       Number of particles simultaneously released at each ', &
2445                    'point: ', I5/)
2446490 FORMAT ('       Particle group ',I2,':'/ &
2447            '          Particle radius: ',E10.3, 'm')
2448491 FORMAT ('          Particle inertia is activated'/ &
[1697]2449            '             density_ratio (rho_fluid/rho_particle) =',F6.3/)
[1]2450492 FORMAT ('          Particles are advected only passively (no inertia)'/)
2451493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
2452            '                                         y:',F8.1,' - ',F8.1,' m'/&
2453            '                                         z:',F8.1,' - ',F8.1,' m'/&
2454            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
2455                       ' m  dz = ',F8.1,' m'/)
2456494 FORMAT ('       Output of particle time series in NetCDF format every ', &
2457                    F8.2,' s'/)
2458495 FORMAT ('       Number of particles in total domain: ',I10/)
[1575]2459496 FORMAT ('       Initial vertical particle positions are interpreted ', &
2460                    'as relative to the given topography')
[1]2461500 FORMAT (//' 1D-Model parameters:'/                           &
2462              ' -------------------'//                           &
2463            '    Simulation time:                   ',F8.1,' s'/ &
2464            '    Run-controll output every:         ',F8.1,' s'/ &
2465            '    Vertical profile output every:     ',F8.1,' s'/ &
2466            '    Mixing length calculation:         ',A/         &
2467            '    Dissipation calculation:           ',A/)
2468502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
[667]2469503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2470504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
[1115]2471505 FORMAT ('    Precipitation parameterization via Seifert-Beheng-Scheme')
[1831]2472506 FORMAT ('    Cloud water sedimentation parameterization via Stokes law')
[1115]2473507 FORMAT ('    Turbulence effects on precipitation process')
2474508 FORMAT ('    Ventilation effects on evaporation of rain drops')
2475509 FORMAT ('    Slope limiter used for sedimentation process')
[1551]2476510 FORMAT ('    Droplet density    :   N_c   = ',F6.1,' 1/cm**3')
2477511 FORMAT ('    Sedimentation Courant number:                  '/&
[1697]2478            '                               C_s   =',F4.1,'        ')
[1429]2479512 FORMAT (/' Date:                 ',A8,6X,'Run:       ',A20/      &
2480            ' Time:                 ',A8,6X,'Run-No.:   ',I2.2/     &
2481            ' Run on host:        ',A10,6X,'En-No.:    ',I2.2)
[1791]2482600 FORMAT (/' Nesting informations:'/ &
2483            ' --------------------'/ &
[1797]2484            ' Nesting mode:                     ',A/ &
2485            ' Nesting-datatransfer mode:        ',A// &
[1791]2486            ' Nest id  parent  number   lower left coordinates   name'/ &
2487            ' (*=me)     id    of PEs      x (m)     y (m)' )
2488601 FORMAT (2X,A1,1X,I2.2,6X,I2.2,5X,I5,5X,F8.2,2X,F8.2,5X,A)
[1]2489
2490 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.