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

Last change on this file since 4178 was 4168, checked in by suehring, 5 years ago

Replace get_topography_top_index functions by pre-calculated arrays in order to save computational resources

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