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

Last change on this file since 3524 was 3524, checked in by raasch, 5 years ago

unused variables removed, missing working precision added, missing preprocessor directives added, bugfix concerning allocation of t_surf_wall_v in nopointer case, declaration statements rearranged to avoid compile time errors, mpi_abort arguments replaced to avoid compile errors

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