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

Last change on this file since 3298 was 3298, checked in by kanani, 6 years ago

Merge chemistry branch at r3297 to trunk

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