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

Last change on this file since 3379 was 3355, checked in by knoop, 6 years ago

removed calc_radiation.f90 and related cloud_top_radiation namelist parameter (functionality replaced by radiation model)

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