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

Last change on this file since 3835 was 3761, checked in by raasch, 6 years ago

unused variables removed, OpenACC directives re-formatted, statements added to avoid compiler warnings

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