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

Last change on this file since 2253 was 2233, checked in by suehring, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 90.2 KB
RevLine 
[1682]1!> @file header.f90
[2000]2!------------------------------------------------------------------------------!
[1036]3! This file is part of PALM.
4!
[2000]5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
[1036]9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
[2101]17! Copyright 1997-2017 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[254]20! Current revisions:
[1]21! -----------------
[1932]22!
[2233]23!
[1485]24! Former revisions:
25! -----------------
26! $Id: header.f90 2233 2017-05-30 18:08:54Z knoop $
27!
[2233]28! 2232 2017-05-30 17:47:52Z suehring
29! Adjustments to new topography and surface concept
30! Generic tunnel setup added
31!
[2201]32! 2200 2017-04-11 11:37:51Z suehring
33! monotonic_adjustment removed
34!
[2119]35! 2118 2017-01-17 16:38:49Z raasch
36! OpenACC relatec code removed
37!
[2074]38! 2073 2016-11-30 14:34:05Z raasch
39! small bugfix concerning output of scalar profiles
40!
[2051]41! 2050 2016-11-08 15:00:55Z gronemeier
42! Implement turbulent outflow condition
43!
[2038]44! 2037 2016-10-26 11:15:40Z knoop
45! Anelastic approximation implemented
46!
[2001]47! 2000 2016-08-20 18:09:15Z knoop
48! Forced header and separation lines into 80 columns
49!
[1993]50! 1992 2016-08-12 15:14:59Z suehring
51! Adapted for top_scalarflux
52!
[1961]53! 1960 2016-07-12 16:34:24Z suehring
54! Treat humidity and passive scalar separately.
55! Modify misleading information concerning humidity.
56! Bugfix, change unit for humidity flux.
57!
[1958]58! 1957 2016-07-07 10:43:48Z suehring
59! flight module added
60!
[1932]61! 1931 2016-06-10 12:06:59Z suehring
62! Rename multigrid into multigrid_noopt
63!
[1903]64! 1902 2016-05-09 11:18:56Z suehring
65! Write information about masking_method only for multigrid solver
66!
[1851]67! 1849 2016-04-08 11:33:18Z hoffmann
68! Adapted for modularization of microphysics
69!
[1834]70! 1833 2016-04-07 14:23:03Z raasch
71! spectrum renamed spectra_mod, output of spectra related quantities moved to
72! spectra_mod
73!
[1832]74! 1831 2016-04-07 13:15:51Z hoffmann
75! turbulence renamed collision_turbulence,
76! drizzle renamed cloud_water_sedimentation
77!
[1827]78! 1826 2016-04-07 12:01:39Z maronga
79! Moved radiation model header output to the respective module.
80! Moved canopy model header output to the respective module.
81!
[1823]82! 1822 2016-04-07 07:49:42Z hoffmann
83! Tails removed. icloud_scheme replaced by microphysics_*
84!
[1818]85! 1817 2016-04-06 15:44:20Z maronga
86! Moved land_surface_model header output to the respective module.
87!
[1809]88! 1808 2016-04-05 19:44:00Z raasch
89! routine local_flush replaced by FORTRAN statement
90!
[1798]91! 1797 2016-03-21 16:50:28Z raasch
92! output of nesting datatransfer mode
93!
[1792]94! 1791 2016-03-11 10:41:25Z raasch
95! output of nesting informations of all domains
96!
[1789]97! 1788 2016-03-10 11:01:04Z maronga
98! Parameter dewfall removed
99!
[1787]100! 1786 2016-03-08 05:49:27Z raasch
101! cpp-direktives for spectra removed
102!
[1784]103! 1783 2016-03-06 18:36:17Z raasch
104! netcdf module and variable names changed, output of netcdf_deflate
105!
[1765]106! 1764 2016-02-28 12:45:19Z raasch
107! output of nesting informations
108!
[1698]109! 1697 2015-10-28 17:14:10Z raasch
110! small E- and F-FORMAT changes to avoid informative compiler messages about
111! insufficient field width
112!
[1692]113! 1691 2015-10-26 16:17:44Z maronga
114! Renamed prandtl_layer to constant_flux_layer, renames rif_min/rif_max to
115! zeta_min/zeta_max.
116!
[1683]117! 1682 2015-10-07 23:56:08Z knoop
118! Code annotations made doxygen readable
119!
[1676]120! 1675 2015-10-02 08:28:59Z gronemeier
121! Bugfix: Definition of topography grid levels
122!
[1662]123! 1660 2015-09-21 08:15:16Z gronemeier
124! Bugfix: Definition of building/street canyon height if vertical grid stretching
125!         starts below the maximum topography height.
126!
[1591]127! 1590 2015-05-08 13:56:27Z maronga
128! Bugfix: Added TRIM statements for character strings for LSM and radiation code
129!
[1586]130! 1585 2015-04-30 07:05:52Z maronga
131! Further output for radiation model(s).
132!
[1576]133! 1575 2015-03-27 09:56:27Z raasch
134! adjustments for psolver-queries, output of seed_follows_topography
135!
[1561]136! 1560 2015-03-06 10:48:54Z keck
137! output for recycling y shift
138!
[1558]139! 1557 2015-03-05 16:43:04Z suehring
140! output for monotonic limiter
141!
[1552]142! 1551 2015-03-03 14:18:16Z maronga
143! Added informal output for land surface model and radiation model. Removed typo.
144!
[1497]145! 1496 2014-12-02 17:25:50Z maronga
146! Renamed: "radiation -> "cloud_top_radiation"
147!
[1485]148! 1484 2014-10-21 10:53:05Z kanani
[1484]149! Changes due to new module structure of the plant canopy model:
150!   module plant_canopy_model_mod and output for new canopy model parameters
151!   (alpha_lad, beta_lad, lai_beta,...) added,
152!   drag_coefficient, leaf_surface_concentration and scalar_exchange_coefficient
153!   renamed to canopy_drag_coeff, leaf_surface_conc and leaf_scalar_exch_coeff,
154!   learde renamed leaf_area_density.
155! Bugfix: DO-WHILE-loop for lad header information additionally restricted
156! by maximum number of gradient levels (currently 10)
[1483]157!
158! 1482 2014-10-18 12:34:45Z raasch
159! information about calculated or predefined virtual processor topology adjusted
160!
[1469]161! 1468 2014-09-24 14:06:57Z maronga
162! Adapted for use on up to 6-digit processor cores
163!
[1430]164! 1429 2014-07-15 12:53:45Z knoop
165! header exended to provide ensemble_member_nr if specified
166!
[1377]167! 1376 2014-04-26 11:21:22Z boeske
168! Correction of typos
169!
[1366]170! 1365 2014-04-22 15:03:56Z boeske
171! New section 'Large scale forcing and nudging':
172! output of large scale forcing and nudging information,
173! new section for initial profiles created
174!
[1360]175! 1359 2014-04-11 17:15:14Z hoffmann
176! dt_sort_particles removed
177!
[1354]178! 1353 2014-04-08 15:21:23Z heinze
179! REAL constants provided with KIND-attribute
180!
[1329]181! 1327 2014-03-21 11:00:16Z raasch
182! parts concerning iso2d and avs output removed,
183! -netcdf output queries
184!
[1325]185! 1324 2014-03-21 09:13:16Z suehring
186! Bugfix: module spectrum added
187!
[1323]188! 1322 2014-03-20 16:38:49Z raasch
189! REAL functions provided with KIND-attribute,
190! some REAL constants defined as wp-kind
191!
[1321]192! 1320 2014-03-20 08:40:49Z raasch
[1320]193! ONLY-attribute added to USE-statements,
194! kind-parameters added to all INTEGER and REAL declaration statements,
195! kinds are defined in new module kinds,
196! revision history before 2012 removed,
197! comment fields (!:) to be used for variable explanations added to
198! all variable declaration statements
[1321]199!
[1309]200! 1308 2014-03-13 14:58:42Z fricke
201! output of the fixed number of output time levels
202! output_format adjusted for masked data if netcdf_data_format > 5
203!
[1300]204! 1299 2014-03-06 13:15:21Z heinze
205! output for using large_scale subsidence in combination
206! with large_scale_forcing
207! reformatting, more detailed explanations
208!
[1242]209! 1241 2013-10-30 11:36:58Z heinze
210! output for nudging + large scale forcing from external file
211!
[1217]212! 1216 2013-08-26 09:31:42Z raasch
213! output for transpose_compute_overlap
214!
[1213]215! 1212 2013-08-15 08:46:27Z raasch
216! output for poisfft_hybrid removed
217!
[1182]218! 1179 2013-06-14 05:57:58Z raasch
219! output of reference_state, use_reference renamed use_single_reference_value
220!
[1160]221! 1159 2013-05-21 11:58:22Z fricke
222! +use_cmax
223!
[1116]224! 1115 2013-03-26 18:16:16Z hoffmann
225! descriptions for Seifert-Beheng-cloud-physics-scheme added
226!
[1112]227! 1111 2013-03-08 23:54:10Z raasch
228! output of accelerator board information
229! ibc_p_b = 2 removed
230!
[1109]231! 1108 2013-03-05 07:03:32Z raasch
232! bugfix for r1106
233!
[1107]234! 1106 2013-03-04 05:31:38Z raasch
235! some format changes for coupled runs
236!
[1093]237! 1092 2013-02-02 11:24:22Z raasch
238! unused variables removed
239!
[1037]240! 1036 2012-10-22 13:43:42Z raasch
241! code put under GPL (PALM 3.9)
242!
[1035]243! 1031 2012-10-19 14:35:30Z raasch
244! output of netCDF data format modified
245!
[1017]246! 1015 2012-09-27 09:23:24Z raasch
[1365]247! output of Adjustment of mixing length to the Prandtl mixing length at first
[1017]248! grid point above ground removed
249!
[1004]250! 1003 2012-09-14 14:35:53Z raasch
251! output of information about equal/unequal subdomain size removed
252!
[1002]253! 1001 2012-09-13 14:08:46Z raasch
254! all actions concerning leapfrog- and upstream-spline-scheme removed
255!
[979]256! 978 2012-08-09 08:28:32Z fricke
257! -km_damp_max, outflow_damping_width
258! +pt_damping_factor, pt_damping_width
259! +z0h
260!
[965]261! 964 2012-07-26 09:14:24Z raasch
262! output of profil-related quantities removed
263!
[941]264! 940 2012-07-09 14:31:00Z raasch
265! Output in case of simulations for pure neutral stratification (no pt-equation
266! solved)
267!
[928]268! 927 2012-06-06 19:15:04Z raasch
269! output of masking_method for mg-solver
270!
[869]271! 868 2012-03-28 12:21:07Z raasch
272! translation velocity in Galilean transformation changed to 0.6 * ug
273!
[834]274! 833 2012-02-22 08:55:55Z maronga
275! Adjusted format for leaf area density
276!
[829]277! 828 2012-02-21 12:00:36Z raasch
278! output of dissipation_classes + radius_classes
279!
[826]280! 825 2012-02-19 03:03:44Z raasch
281! Output of cloud physics parameters/quantities complemented and restructured
282!
[1]283! Revision 1.1  1997/08/11 06:17:20  raasch
284! Initial revision
285!
286!
287! Description:
288! ------------
[1764]289!> Writing a header with all important information about the current run.
[1682]290!> This subroutine is called three times, two times at the beginning
291!> (writing information on files RUN_CONTROL and HEADER) and one time at the
292!> end of the run, then writing additional information about CPU-usage on file
293!> header.
[411]294!-----------------------------------------------------------------------------!
[1682]295 SUBROUTINE header
296 
[1]297
[1320]298    USE arrays_3d,                                                             &
[2232]299        ONLY:  pt_init, q_init, s_init, sa_init, ug, vg, w_subs, zu, zw
[1320]300       
[1]301    USE control_parameters
[1320]302       
303    USE cloud_parameters,                                                      &
[1849]304        ONLY:  cp, l_v, r_d
305
[1320]306    USE cpulog,                                                                &
307        ONLY:  log_point_s
308       
309    USE dvrp_variables,                                                        &
310        ONLY:  use_seperate_pe_for_dvrp_output
311       
[1957]312    USE flight_mod,                                                            &
313        ONLY:  flight_header
314       
[1320]315    USE grid_variables,                                                        &
316        ONLY:  dx, dy
317       
318    USE indices,                                                               &
319        ONLY:  mg_loc_ind, nnx, nny, nnz, nx, ny, nxl_mg, nxr_mg, nyn_mg,      &
320               nys_mg, nzt, nzt_mg
321       
322    USE kinds
[1817]323 
[1551]324    USE land_surface_model_mod,                                                &
[2232]325        ONLY: lsm_header
[1849]326
327    USE microphysics_mod,                                                      &
328        ONLY:  cloud_water_sedimentation, collision_turbulence,                &
329               c_sedimentation, limiter_sedimentation, nc_const,               &
330               ventilation_effect
331
[1320]332    USE model_1d,                                                              &
333        ONLY:  damp_level_ind_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
334       
[1783]335    USE netcdf_interface,                                                      &
336        ONLY:  netcdf_data_format, netcdf_data_format_string, netcdf_deflate
337
[1320]338    USE particle_attributes,                                                   &
339        ONLY:  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, collision_kernel,     &
[1831]340               curvature_solution_effects,                                     &
[1320]341               density_ratio, dissipation_classes, dt_min_part, dt_prel,       &
[1359]342               dt_write_particle_data, end_time_prel,                          &
[1822]343               number_of_particle_groups, particle_advection,                  &
344               particle_advection_start,                                       &
[1320]345               particles_per_point, pdx, pdy, pdz,  psb, psl, psn, psr, pss,   &
346               pst, radius, radius_classes, random_start_position,             &
[1575]347               seed_follows_topography,                                        &
[1822]348               total_number_of_particles, use_sgs_for_particles,               &
[1320]349               vertical_particle_advection, write_particle_statistics
350       
[1]351    USE pegrid
[1484]352
353    USE plant_canopy_model_mod,                                                &
[1826]354        ONLY:  pcm_header, plant_canopy
[1551]355
[1791]356    USE pmc_handle_communicator,                                               &
357        ONLY:  pmc_get_model_info
358
[1764]359    USE pmc_interface,                                                         &
[1797]360        ONLY:  nested_run, nesting_datatransfer_mode, nesting_mode
[1764]361
[1551]362    USE radiation_model_mod,                                                   &
[1826]363        ONLY:  radiation, radiation_header
[1324]364   
[1833]365    USE spectra_mod,                                                           &
366        ONLY:  calculate_spectra, spectra_header
[1]367
[2232]368    USE surface_mod,                                                           &
369        ONLY:  surf_def_h
370
[1]371    IMPLICIT NONE
372
[1682]373    CHARACTER (LEN=1)  ::  prec                !<
[1320]374   
[1682]375    CHARACTER (LEN=2)  ::  do2d_mode           !<
[1320]376   
[1682]377    CHARACTER (LEN=5)  ::  section_chr         !<
[1320]378   
[1682]379    CHARACTER (LEN=10) ::  coor_chr            !<
380    CHARACTER (LEN=10) ::  host_chr            !<
[1320]381   
[1682]382    CHARACTER (LEN=16) ::  begin_chr           !<
[1320]383   
[1682]384    CHARACTER (LEN=26) ::  ver_rev             !<
[1791]385
386    CHARACTER (LEN=32) ::  cpl_name            !<
[1320]387   
[1682]388    CHARACTER (LEN=40) ::  output_format       !<
[1320]389   
[1682]390    CHARACTER (LEN=70) ::  char1               !<
391    CHARACTER (LEN=70) ::  char2               !<
392    CHARACTER (LEN=70) ::  dopr_chr            !<
393    CHARACTER (LEN=70) ::  do2d_xy             !<
394    CHARACTER (LEN=70) ::  do2d_xz             !<
395    CHARACTER (LEN=70) ::  do2d_yz             !<
396    CHARACTER (LEN=70) ::  do3d_chr            !<
397    CHARACTER (LEN=70) ::  domask_chr          !<
398    CHARACTER (LEN=70) ::  run_classification  !<
[1320]399   
[1826]400    CHARACTER (LEN=85) ::  r_upper             !<
401    CHARACTER (LEN=85) ::  r_lower             !<
[1320]402   
[1682]403    CHARACTER (LEN=86) ::  coordinates         !<
404    CHARACTER (LEN=86) ::  gradients           !<
405    CHARACTER (LEN=86) ::  slices              !<
406    CHARACTER (LEN=86) ::  temperatures        !<
407    CHARACTER (LEN=86) ::  ugcomponent         !<
408    CHARACTER (LEN=86) ::  vgcomponent         !<
[1]409
[1682]410    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)  !<
[410]411
[1791]412    INTEGER(iwp) ::  av             !<
413    INTEGER(iwp) ::  bh             !<
414    INTEGER(iwp) ::  blx            !<
415    INTEGER(iwp) ::  bly            !<
416    INTEGER(iwp) ::  bxl            !<
417    INTEGER(iwp) ::  bxr            !<
418    INTEGER(iwp) ::  byn            !<
419    INTEGER(iwp) ::  bys            !<
420    INTEGER(iwp) ::  ch             !<
421    INTEGER(iwp) ::  count          !<
422    INTEGER(iwp) ::  cpl_parent_id  !<
423    INTEGER(iwp) ::  cwx            !<
424    INTEGER(iwp) ::  cwy            !<
425    INTEGER(iwp) ::  cxl            !<
426    INTEGER(iwp) ::  cxr            !<
427    INTEGER(iwp) ::  cyn            !<
428    INTEGER(iwp) ::  cys            !<
429    INTEGER(iwp) ::  dim            !<
430    INTEGER(iwp) ::  i              !<
431    INTEGER(iwp) ::  io             !<
432    INTEGER(iwp) ::  j              !<
433    INTEGER(iwp) ::  k              !<
434    INTEGER(iwp) ::  l              !<
435    INTEGER(iwp) ::  ll             !<
436    INTEGER(iwp) ::  mpi_type       !<
437    INTEGER(iwp) ::  my_cpl_id      !<
438    INTEGER(iwp) ::  n              !<
439    INTEGER(iwp) ::  ncpl           !<
440    INTEGER(iwp) ::  npe_total      !<
[1320]441   
[1826]442
[1682]443    REAL(wp) ::  cpuseconds_per_simulated_second  !<
[1791]444    REAL(wp) ::  lower_left_coord_x               !< x-coordinate of nest domain
445    REAL(wp) ::  lower_left_coord_y               !< y-coordinate of nest domain
[1]446
447!
448!-- Open the output file. At the end of the simulation, output is directed
449!-- to unit 19.
450    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
451         .NOT. simulated_time_at_begin /= simulated_time )  THEN
452       io = 15   !  header output on file RUN_CONTROL
453    ELSE
454       io = 19   !  header output on file HEADER
455    ENDIF
456    CALL check_open( io )
457
458!
459!-- At the end of the run, output file (HEADER) will be rewritten with
[1551]460!-- new information
[1]461    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
462
463!
464!-- Determine kind of model run
465    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
[1764]466       run_classification = 'restart run'
[328]467    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[1764]468       run_classification = 'run with cyclic fill of 3D - prerun data'
[147]469    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
[1764]470       run_classification = 'run without 1D - prerun'
[197]471    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
[1764]472       run_classification = 'run with 1D - prerun'
[197]473    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
[1764]474       run_classification = 'run initialized by user'
[1]475    ELSE
[254]476       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
477       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
[1]478    ENDIF
[1764]479    IF ( nested_run )  run_classification = 'nested ' // run_classification
[97]480    IF ( ocean )  THEN
481       run_classification = 'ocean - ' // run_classification
482    ELSE
483       run_classification = 'atmosphere - ' // run_classification
484    ENDIF
[1]485
486!
487!-- Run-identification, date, time, host
488    host_chr = host(1:10)
[75]489    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
[102]490    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
[291]491    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
492#if defined( __mpi2 )
493       mpi_type = 2
494#else
495       mpi_type = 1
496#endif
497       WRITE ( io, 101 )  mpi_type, coupling_mode
498    ENDIF
[1108]499#if defined( __parallel )
[1353]500    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]501       IF ( coupling_start_time > simulated_time_at_begin )  THEN
502          WRITE ( io, 109 )
503       ELSE
504          WRITE ( io, 114 )
505       ENDIF
506    ENDIF
[1108]507#endif
[1429]508    IF ( ensemble_member_nr /= 0 )  THEN
509       WRITE ( io, 512 )  run_date, run_identifier, run_time, runnr,           &
510                       ADJUSTR( host_chr ), ensemble_member_nr
511    ELSE
512       WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr,           &
[102]513                       ADJUSTR( host_chr )
[1429]514    ENDIF
[1]515#if defined( __parallel )
[1482]516    IF ( npex == -1  .AND.  npey == -1 )  THEN
[1]517       char1 = 'calculated'
518    ELSE
519       char1 = 'predefined'
520    ENDIF
521    IF ( threads_per_task == 1 )  THEN
[102]522       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
[1]523    ELSE
[102]524       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
[1]525                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
526    ENDIF
527    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
528           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
529         npex == -1  .AND.  pdims(2) == 1 )                      &
530    THEN
[102]531       WRITE ( io, 106 )
[1]532    ELSEIF ( pdims(2) == 1 )  THEN
[102]533       WRITE ( io, 107 )  'x'
[1]534    ELSEIF ( pdims(1) == 1 )  THEN
[102]535       WRITE ( io, 107 )  'y'
[1]536    ENDIF
[102]537    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
[759]538    IF ( numprocs /= maximum_parallel_io_streams )  THEN
539       WRITE ( io, 108 )  maximum_parallel_io_streams
540    ENDIF
[1]541#endif
[1764]542
543!
544!-- Nesting informations
545    IF ( nested_run )  THEN
[1791]546
[1797]547       WRITE ( io, 600 )  TRIM( nesting_mode ),                                &
548                          TRIM( nesting_datatransfer_mode )
[1791]549       CALL pmc_get_model_info( ncpl = ncpl, cpl_id = my_cpl_id )
550
551       DO  n = 1, ncpl
552          CALL pmc_get_model_info( request_for_cpl_id = n, cpl_name = cpl_name,&
553                                   cpl_parent_id = cpl_parent_id,              &
554                                   lower_left_x = lower_left_coord_x,          &
555                                   lower_left_y = lower_left_coord_y,          &
556                                   npe_total = npe_total )
557          IF ( n == my_cpl_id )  THEN
558             char1 = '*'
559          ELSE
560             char1 = ' '
561          ENDIF
562          WRITE ( io, 601 )  TRIM( char1 ), n, cpl_parent_id, npe_total,       &
563                             lower_left_coord_x, lower_left_coord_y,           &
564                             TRIM( cpl_name )
565       ENDDO
[1764]566    ENDIF
[1]567    WRITE ( io, 99 )
568
569!
570!-- Numerical schemes
571    WRITE ( io, 110 )
[2037]572    WRITE ( io, 121 )  TRIM( approximation )
[1]573    IF ( psolver(1:7) == 'poisfft' )  THEN
574       WRITE ( io, 111 )  TRIM( fft_method )
[1216]575       IF ( transpose_compute_overlap )  WRITE( io, 115 )
[1]576    ELSEIF ( psolver == 'sor' )  THEN
577       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
[1575]578    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
579       WRITE ( io, 135 )  TRIM(psolver), cycle_mg, maximum_grid_level, ngsrb
[1]580       IF ( mg_cycles == -1 )  THEN
581          WRITE ( io, 140 )  residual_limit
582       ELSE
583          WRITE ( io, 141 )  mg_cycles
584       ENDIF
585       IF ( mg_switch_to_pe0_level == 0 )  THEN
586          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
587                             nzt_mg(1)
[197]588       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
[1]589          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
590                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
591                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
592                             nzt_mg(mg_switch_to_pe0_level),    &
593                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
594                             nzt_mg(1)
595       ENDIF
[1931]596       IF ( psolver == 'multigrid_noopt' .AND. masking_method )  WRITE ( io, 144 )
[1]597    ENDIF
598    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
599    THEN
600       WRITE ( io, 142 )
601    ENDIF
602
603    IF ( momentum_advec == 'pw-scheme' )  THEN
604       WRITE ( io, 113 )
[1299]605    ELSEIF (momentum_advec == 'ws-scheme' )  THEN
[667]606       WRITE ( io, 503 )
[1]607    ENDIF
608    IF ( scalar_advec == 'pw-scheme' )  THEN
609       WRITE ( io, 116 )
[667]610    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
611       WRITE ( io, 504 )
[1]612    ELSE
613       WRITE ( io, 118 )
614    ENDIF
[63]615
616    WRITE ( io, 139 )  TRIM( loop_optimization )
617
[1]618    IF ( galilei_transformation )  THEN
619       IF ( use_ug_for_galilei_tr )  THEN
[868]620          char1 = '0.6 * geostrophic wind'
[1]621       ELSE
622          char1 = 'mean wind in model domain'
623       ENDIF
624       IF ( simulated_time_at_begin == simulated_time )  THEN
625          char2 = 'at the start of the run'
626       ELSE
627          char2 = 'at the end of the run'
628       ENDIF
[1353]629       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ),                        &
630                          advected_distance_x/1000.0_wp,                       &
631                          advected_distance_y/1000.0_wp
[1]632    ENDIF
[1001]633    WRITE ( io, 122 )  timestep_scheme
[87]634    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
[1353]635    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
[108]636       IF ( .NOT. ocean )  THEN
637          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
638               rayleigh_damping_factor
639       ELSE
640          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
641               rayleigh_damping_factor
642       ENDIF
[1]643    ENDIF
[940]644    IF ( neutral )  WRITE ( io, 131 )  pt_surface
[75]645    IF ( humidity )  THEN
[1]646       IF ( .NOT. cloud_physics )  THEN
647          WRITE ( io, 129 )
648       ELSE
649          WRITE ( io, 130 )
650       ENDIF
651    ENDIF
652    IF ( passive_scalar )  WRITE ( io, 134 )
[240]653    IF ( conserve_volume_flow )  THEN
[241]654       WRITE ( io, 150 )  conserve_volume_flow_mode
655       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
656          WRITE ( io, 151 )  u_bulk, v_bulk
657       ENDIF
[240]658    ELSEIF ( dp_external )  THEN
659       IF ( dp_smooth )  THEN
[241]660          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
[240]661       ELSE
[241]662          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
[240]663       ENDIF
664    ENDIF
[1]665    WRITE ( io, 99 )
666
667!
[1551]668!-- Runtime and timestep information
[1]669    WRITE ( io, 200 )
670    IF ( .NOT. dt_fixed )  THEN
671       WRITE ( io, 201 )  dt_max, cfl_factor
672    ELSE
673       WRITE ( io, 202 )  dt
674    ENDIF
675    WRITE ( io, 203 )  simulated_time_at_begin, end_time
676
[1322]677    IF ( time_restart /= 9999999.9_wp  .AND. &
[1]678         simulated_time_at_begin == simulated_time )  THEN
[1322]679       IF ( dt_restart == 9999999.9_wp )  THEN
[1]680          WRITE ( io, 204 )  ' Restart at:       ',time_restart
681       ELSE
682          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
683       ENDIF
684    ENDIF
685
686    IF ( simulated_time_at_begin /= simulated_time )  THEN
687       i = MAX ( log_point_s(10)%counts, 1 )
[1353]688       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0_wp )  THEN
689          cpuseconds_per_simulated_second = 0.0_wp
[1]690       ELSE
691          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
692                                            ( simulated_time -    &
693                                              simulated_time_at_begin )
694       ENDIF
[1322]695       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum,      &
696                          log_point_s(10)%sum / REAL( i, KIND=wp ), &
[1]697                          cpuseconds_per_simulated_second
[1322]698       IF ( time_restart /= 9999999.9_wp  .AND.  time_restart < end_time )  THEN
699          IF ( dt_restart == 9999999.9_wp )  THEN
[1106]700             WRITE ( io, 204 )  ' Next restart at:     ',time_restart
[1]701          ELSE
[1106]702             WRITE ( io, 205 )  ' Next restart at:     ',time_restart, dt_restart
[1]703          ENDIF
704       ENDIF
705    ENDIF
706
[1324]707
[1]708!
[291]709!-- Start time for coupled runs, if independent precursor runs for atmosphere
[1106]710!-- and ocean are used or have been used. In this case, coupling_start_time
711!-- defines the time when the coupling is switched on.
[1353]712    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]713       WRITE ( io, 207 )  coupling_start_time
[291]714    ENDIF
715
716!
[1]717!-- Computational grid
[94]718    IF ( .NOT. ocean )  THEN
719       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
720       IF ( dz_stretch_level_index < nzt+1 )  THEN
721          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
722                             dz_stretch_factor, dz_max
723       ENDIF
724    ELSE
725       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
726       IF ( dz_stretch_level_index > 0 )  THEN
727          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
728                             dz_stretch_factor, dz_max
729       ENDIF
[1]730    ENDIF
731    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
732                       MIN( nnz+2, nzt+2 )
733    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
734
735!
[1365]736!-- Large scale forcing and nudging
737    WRITE ( io, 160 )
738    IF ( large_scale_forcing )  THEN
739       WRITE ( io, 162 )
740       WRITE ( io, 163 )
741
742       IF ( large_scale_subsidence )  THEN
743          IF ( .NOT. use_subsidence_tendencies )  THEN
744             WRITE ( io, 164 )
745          ELSE
746             WRITE ( io, 165 )
747          ENDIF
748       ENDIF
749
750       IF ( bc_pt_b == 'dirichlet' )  THEN
751          WRITE ( io, 180 )
752       ELSEIF ( bc_pt_b == 'neumann' )  THEN
753          WRITE ( io, 181 )
754       ENDIF
755
756       IF ( bc_q_b == 'dirichlet' )  THEN
757          WRITE ( io, 182 )
758       ELSEIF ( bc_q_b == 'neumann' )  THEN
759          WRITE ( io, 183 )
760       ENDIF
761
762       WRITE ( io, 167 )
763       IF ( nudging )  THEN
764          WRITE ( io, 170 )
765       ENDIF
766    ELSE
767       WRITE ( io, 161 )
768       WRITE ( io, 171 )
769    ENDIF
770    IF ( large_scale_subsidence )  THEN
771       WRITE ( io, 168 )
772       WRITE ( io, 169 )
773    ENDIF
774
775!
776!-- Profile for the large scale vertial velocity
777!-- Building output strings, starting with surface value
778    IF ( large_scale_subsidence )  THEN
779       temperatures = '   0.0'
780       gradients = '------'
781       slices = '     0'
782       coordinates = '   0.0'
783       i = 1
784       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
785
786          WRITE (coor_chr,'(E10.2,7X)')  &
787                                w_subs(subs_vertical_gradient_level_i(i))
788          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
789
790          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
791          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
792
793          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
794          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
795
796          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
797          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
798
799          IF ( i == 10 )  THEN
800             EXIT
801          ELSE
802             i = i + 1
803          ENDIF
804
805       ENDDO
806
807 
808       IF ( .NOT. large_scale_forcing )  THEN
809          WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
810                             TRIM( gradients ), TRIM( slices )
811       ENDIF
812
813
814    ENDIF
815
816!-- Profile of the geostrophic wind (component ug)
817!-- Building output strings
818    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
819    gradients = '------'
820    slices = '     0'
821    coordinates = '   0.0'
822    i = 1
823    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
824     
825       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
826       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
827
828       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
829       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
830
831       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
832       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
833
834       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
835       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
836
837       IF ( i == 10 )  THEN
838          EXIT
839       ELSE
840          i = i + 1
841       ENDIF
842
843    ENDDO
844
845    IF ( .NOT. large_scale_forcing )  THEN
846       WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
847                          TRIM( gradients ), TRIM( slices )
848    ENDIF
849
850!-- Profile of the geostrophic wind (component vg)
851!-- Building output strings
852    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
853    gradients = '------'
854    slices = '     0'
855    coordinates = '   0.0'
856    i = 1
857    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
858
859       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
860       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
861
862       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
863       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
864
865       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
866       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
867
868       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
869       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
870
871       IF ( i == 10 )  THEN
872          EXIT
873       ELSE
874          i = i + 1
875       ENDIF
876 
877    ENDDO
878
879    IF ( .NOT. large_scale_forcing )  THEN
880       WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
881                          TRIM( gradients ), TRIM( slices )
882    ENDIF
883
884!
[1]885!-- Topography
886    WRITE ( io, 270 )  topography
887    SELECT CASE ( TRIM( topography ) )
888
889       CASE ( 'flat' )
890          ! no actions necessary
891
892       CASE ( 'single_building' )
893          blx = INT( building_length_x / dx )
894          bly = INT( building_length_y / dy )
[1675]895          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
896          IF ( ABS( zw(bh  ) - building_height ) == &
897               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
[1]898
[1322]899          IF ( building_wall_left == 9999999.9_wp )  THEN
[1]900             building_wall_left = ( nx + 1 - blx ) / 2 * dx
901          ENDIF
[1353]902          bxl = INT ( building_wall_left / dx + 0.5_wp )
[1]903          bxr = bxl + blx
904
[1322]905          IF ( building_wall_south == 9999999.9_wp )  THEN
[1]906             building_wall_south = ( ny + 1 - bly ) / 2 * dy
907          ENDIF
[1353]908          bys = INT ( building_wall_south / dy + 0.5_wp )
[1]909          byn = bys + bly
910
911          WRITE ( io, 271 )  building_length_x, building_length_y, &
912                             building_height, bxl, bxr, bys, byn
913
[240]914       CASE ( 'single_street_canyon' )
[1675]915          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
916          IF ( ABS( zw(ch  ) - canyon_height ) == &
917               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
[1322]918          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[240]919!
920!--          Street canyon in y direction
921             cwx = NINT( canyon_width_x / dx )
[1322]922             IF ( canyon_wall_left == 9999999.9_wp )  THEN
[240]923                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
924             ENDIF
925             cxl = NINT( canyon_wall_left / dx )
926             cxr = cxl + cwx
927             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
928
[1322]929          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[240]930!
931!--          Street canyon in x direction
932             cwy = NINT( canyon_width_y / dy )
[1322]933             IF ( canyon_wall_south == 9999999.9_wp )  THEN
[240]934                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
935             ENDIF
936             cys = NINT( canyon_wall_south / dy )
937             cyn = cys + cwy
938             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
939          ENDIF
940
[2232]941       CASE ( 'tunnel' )
942          IF ( tunnel_width_x /= 9999999.9_wp )  THEN
943!
944!--          Tunnel axis in y direction
945             IF ( tunnel_length == 9999999.9_wp  .OR.                          &
946                  tunnel_length >= ( nx + 1 ) * dx )  THEN
947                WRITE ( io, 273 )  'y', tunnel_height, tunnel_wall_depth,      &
948                                        tunnel_width_x
949             ELSE
950                WRITE ( io, 274 )  'y', tunnel_height, tunnel_wall_depth,      &
951                                        tunnel_width_x, tunnel_length
952             ENDIF
953
954          ELSEIF ( tunnel_width_y /= 9999999.9_wp )  THEN
955!
956!--          Tunnel axis in x direction
957             IF ( tunnel_length == 9999999.9_wp  .OR.                          &
958                  tunnel_length >= ( ny + 1 ) * dy )  THEN
959                WRITE ( io, 273 )  'x', tunnel_height, tunnel_wall_depth,      &
960                                        tunnel_width_y
961             ELSE
962                WRITE ( io, 274 )  'x', tunnel_height, tunnel_wall_depth,      &
963                                        tunnel_width_y, tunnel_length
964             ENDIF
965          ENDIF
966
[1]967    END SELECT
968
[256]969    IF ( TRIM( topography ) /= 'flat' )  THEN
970       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
971          IF ( TRIM( topography ) == 'single_building' .OR.  &
972               TRIM( topography ) == 'single_street_canyon' )  THEN
973             WRITE ( io, 278 )
974          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
975             WRITE ( io, 279 )
976          ENDIF
977       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
978          WRITE ( io, 278 )
979       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
980          WRITE ( io, 279 )
981       ENDIF
982    ENDIF
983
[1826]984    IF ( plant_canopy )  CALL pcm_header ( io )
[138]985
[1817]986    IF ( land_surface )  CALL lsm_header ( io )
[1484]987
[1826]988    IF ( radiation )  CALL radiation_header ( io )
[1551]989
990!
[1]991!-- Boundary conditions
992    IF ( ibc_p_b == 0 )  THEN
[1826]993       r_lower = 'p(0)     = 0      |'
[1]994    ELSEIF ( ibc_p_b == 1 )  THEN
[1826]995       r_lower = 'p(0)     = p(1)   |'
[1]996    ENDIF
997    IF ( ibc_p_t == 0 )  THEN
[1826]998       r_upper  = 'p(nzt+1) = 0      |'
[1]999    ELSE
[1826]1000       r_upper  = 'p(nzt+1) = p(nzt) |'
[1]1001    ENDIF
1002
1003    IF ( ibc_uv_b == 0 )  THEN
[1826]1004       r_lower = TRIM( r_lower ) // ' uv(0)     = -uv(1)                |'
[1]1005    ELSE
[1826]1006       r_lower = TRIM( r_lower ) // ' uv(0)     = uv(1)                 |'
[1]1007    ENDIF
[132]1008    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
[1826]1009       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = 0                     |'
[132]1010    ELSEIF ( ibc_uv_t == 0 )  THEN
[1826]1011       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
[1]1012    ELSE
[1826]1013       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = uv(nzt)               |'
[1]1014    ENDIF
1015
1016    IF ( ibc_pt_b == 0 )  THEN
[1551]1017       IF ( land_surface )  THEN
[1826]1018          r_lower = TRIM( r_lower ) // ' pt(0)     = from soil model'
[1551]1019       ELSE
[1826]1020          r_lower = TRIM( r_lower ) // ' pt(0)     = pt_surface'
[1551]1021       ENDIF
[102]1022    ELSEIF ( ibc_pt_b == 1 )  THEN
[1826]1023       r_lower = TRIM( r_lower ) // ' pt(0)     = pt(1)'
[102]1024    ELSEIF ( ibc_pt_b == 2 )  THEN
[1826]1025       r_lower = TRIM( r_lower ) // ' pt(0)     = from coupled model'
[1]1026    ENDIF
1027    IF ( ibc_pt_t == 0 )  THEN
[1826]1028       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt_top'
[19]1029    ELSEIF( ibc_pt_t == 1 )  THEN
[1826]1030       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt)'
[19]1031    ELSEIF( ibc_pt_t == 2 )  THEN
[1826]1032       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
[667]1033
[1]1034    ENDIF
1035
[1826]1036    WRITE ( io, 300 )  r_lower, r_upper
[1]1037
1038    IF ( .NOT. constant_diffusion )  THEN
1039       IF ( ibc_e_b == 1 )  THEN
[1826]1040          r_lower = 'e(0)     = e(1)'
[1]1041       ELSE
[1826]1042          r_lower = 'e(0)     = e(1) = (u*/0.1)**2'
[1]1043       ENDIF
[1826]1044       r_upper = 'e(nzt+1) = e(nzt) = e(nzt-1)'
[1]1045
[1826]1046       WRITE ( io, 301 )  'e', r_lower, r_upper       
[1]1047
1048    ENDIF
1049
[97]1050    IF ( ocean )  THEN
[1826]1051       r_lower = 'sa(0)    = sa(1)'
[97]1052       IF ( ibc_sa_t == 0 )  THEN
[1826]1053          r_upper =  'sa(nzt+1) = sa_surface'
[1]1054       ELSE
[1826]1055          r_upper =  'sa(nzt+1) = sa(nzt)'
[1]1056       ENDIF
[1826]1057       WRITE ( io, 301 ) 'sa', r_lower, r_upper
[97]1058    ENDIF
[1]1059
[97]1060    IF ( humidity )  THEN
1061       IF ( ibc_q_b == 0 )  THEN
[1551]1062          IF ( land_surface )  THEN
[1826]1063             r_lower = 'q(0)     = from soil model'
[1551]1064          ELSE
[1826]1065             r_lower = 'q(0)     = q_surface'
[1551]1066          ENDIF
1067
[97]1068       ELSE
[1992]1069          r_lower = 'q(0)      = q(1)'
[97]1070       ENDIF
1071       IF ( ibc_q_t == 0 )  THEN
[1992]1072          r_upper =  'q(nzt+1) = q_top'
[97]1073       ELSE
[1992]1074          r_upper =  'q(nzt+1) = q(nzt) + dq/dz'
[97]1075       ENDIF
[1826]1076       WRITE ( io, 301 ) 'q', r_lower, r_upper
[97]1077    ENDIF
[1]1078
[97]1079    IF ( passive_scalar )  THEN
[1960]1080       IF ( ibc_s_b == 0 )  THEN
[1992]1081          r_lower = 's(0)      = s_surface'
[97]1082       ELSE
[1992]1083          r_lower = 's(0)      = s(1)'
[97]1084       ENDIF
[1960]1085       IF ( ibc_s_t == 0 )  THEN
[1992]1086          r_upper =  's(nzt+1) = s_top'
1087       ELSEIF ( ibc_s_t == 1 )  THEN
1088          r_upper =  's(nzt+1) = s(nzt)'
1089       ELSEIF ( ibc_s_t == 2 )  THEN
1090          r_upper =  's(nzt+1) = s(nzt) + ds/dz'
[97]1091       ENDIF
[1826]1092       WRITE ( io, 301 ) 's', r_lower, r_upper
[1]1093    ENDIF
1094
1095    IF ( use_surface_fluxes )  THEN
1096       WRITE ( io, 303 )
1097       IF ( constant_heatflux )  THEN
[1299]1098          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[2232]1099             IF ( surf_def_h(0)%ns >= 1 )  WRITE ( io, 306 )  surf_def_h(0)%shf(1)
[1241]1100          ELSE
1101             WRITE ( io, 306 )  surface_heatflux
1102          ENDIF
[1]1103          IF ( random_heatflux )  WRITE ( io, 307 )
1104       ENDIF
[75]1105       IF ( humidity  .AND.  constant_waterflux )  THEN
[1299]1106          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[2232]1107             WRITE ( io, 311 ) surf_def_h(0)%qsws(1)
[1241]1108          ELSE
1109             WRITE ( io, 311 ) surface_waterflux
1110          ENDIF
[1]1111       ENDIF
[1960]1112       IF ( passive_scalar  .AND.  constant_scalarflux )  THEN
1113          WRITE ( io, 313 ) surface_scalarflux
[1]1114       ENDIF
1115    ENDIF
1116
[19]1117    IF ( use_top_fluxes )  THEN
1118       WRITE ( io, 304 )
[102]1119       IF ( coupling_mode == 'uncoupled' )  THEN
[151]1120          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
[102]1121          IF ( constant_top_heatflux )  THEN
1122             WRITE ( io, 306 )  top_heatflux
1123          ENDIF
1124       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1125          WRITE ( io, 316 )
[19]1126       ENDIF
[1992]1127       IF ( ocean  .AND.  constant_top_salinityflux )                          &
[97]1128          WRITE ( io, 309 )  top_salinityflux
[1960]1129       IF ( humidity       )  WRITE ( io, 315 )
[1992]1130       IF ( passive_scalar .AND.  constant_top_scalarflux )                    &
1131          WRITE ( io, 302 ) top_scalarflux
[19]1132    ENDIF
1133
[1691]1134    IF ( constant_flux_layer )  THEN
1135       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length,                     &
1136                          z0h_factor*roughness_length, kappa,                  &
1137                          zeta_min, zeta_max
[1]1138       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
[75]1139       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
[1]1140          WRITE ( io, 312 )
1141       ENDIF
[1960]1142       IF ( passive_scalar  .AND.  .NOT. constant_scalarflux )  THEN
[1]1143          WRITE ( io, 314 )
1144       ENDIF
1145    ELSE
1146       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
[1691]1147          WRITE ( io, 310 )  zeta_min, zeta_max
[1]1148       ENDIF
1149    ENDIF
1150
1151    WRITE ( io, 317 )  bc_lr, bc_ns
[707]1152    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1159]1153       WRITE ( io, 318 )  use_cmax, pt_damping_width, pt_damping_factor       
[151]1154       IF ( turbulent_inflow )  THEN
[1560]1155          IF ( .NOT. recycling_yshift ) THEN
1156             WRITE ( io, 319 )  recycling_width, recycling_plane, &
1157                                inflow_damping_height, inflow_damping_width
1158          ELSE
1159             WRITE ( io, 322 )  recycling_width, recycling_plane, &
1160                                inflow_damping_height, inflow_damping_width
1161          END IF
[151]1162       ENDIF
[2050]1163       IF ( turbulent_outflow )  THEN
1164          WRITE ( io, 323 )  outflow_source_plane, INT(outflow_source_plane/dx)
1165       ENDIF
[1]1166    ENDIF
1167
1168!
[1365]1169!-- Initial Profiles
1170    WRITE ( io, 321 )
1171!
1172!-- Initial wind profiles
1173    IF ( u_profile(1) /= 9999999.9_wp )  WRITE ( io, 427 )
1174
1175!
1176!-- Initial temperature profile
1177!-- Building output strings, starting with surface temperature
1178    WRITE ( temperatures, '(F6.2)' )  pt_surface
1179    gradients = '------'
1180    slices = '     0'
1181    coordinates = '   0.0'
1182    i = 1
1183    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1184
1185       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1186       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1187
1188       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1189       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1190
1191       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1192       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1193
1194       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1195       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1196
1197       IF ( i == 10 )  THEN
1198          EXIT
1199       ELSE
1200          i = i + 1
1201       ENDIF
1202
1203    ENDDO
1204
1205    IF ( .NOT. nudging )  THEN
1206       WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1207                          TRIM( gradients ), TRIM( slices )
1208    ELSE
1209       WRITE ( io, 428 ) 
1210    ENDIF
1211
1212!
1213!-- Initial humidity profile
1214!-- Building output strings, starting with surface humidity
[1960]1215    IF ( humidity )  THEN
[1365]1216       WRITE ( temperatures, '(E8.1)' )  q_surface
1217       gradients = '--------'
1218       slices = '       0'
1219       coordinates = '     0.0'
1220       i = 1
1221       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1222         
1223          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1224          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1225
1226          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1227          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1228         
1229          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1230          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1231         
1232          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1233          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1234
1235          IF ( i == 10 )  THEN
1236             EXIT
1237          ELSE
1238             i = i + 1
1239          ENDIF
1240
1241       ENDDO
1242
[1960]1243       IF ( .NOT. nudging )  THEN
1244          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ),        &
[1365]1245                             TRIM( gradients ), TRIM( slices )
1246       ENDIF
1247    ENDIF
[1960]1248!
1249!-- Initial scalar profile
1250!-- Building output strings, starting with surface humidity
1251    IF ( passive_scalar )  THEN
1252       WRITE ( temperatures, '(E8.1)' )  s_surface
1253       gradients = '--------'
1254       slices = '       0'
1255       coordinates = '     0.0'
1256       i = 1
1257       DO  WHILE ( s_vertical_gradient_level_ind(i) /= -9999 )
1258         
[2073]1259          WRITE (coor_chr,'(E8.1,4X)')  s_init(s_vertical_gradient_level_ind(i))
[1960]1260          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
[1365]1261
[1960]1262          WRITE (coor_chr,'(E8.1,4X)')  s_vertical_gradient(i)
1263          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1264         
1265          WRITE (coor_chr,'(I8,4X)')  s_vertical_gradient_level_ind(i)
1266          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1267         
1268          WRITE (coor_chr,'(F8.1,4X)')  s_vertical_gradient_level(i)
1269          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1270
1271          IF ( i == 10 )  THEN
1272             EXIT
1273          ELSE
1274             i = i + 1
1275          ENDIF
1276
1277       ENDDO
1278
1279       WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ),           &
1280                          TRIM( gradients ), TRIM( slices )
1281    ENDIF   
1282
[1365]1283!
1284!-- Initial salinity profile
1285!-- Building output strings, starting with surface salinity
1286    IF ( ocean )  THEN
1287       WRITE ( temperatures, '(F6.2)' )  sa_surface
1288       gradients = '------'
1289       slices = '     0'
1290       coordinates = '   0.0'
1291       i = 1
1292       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1293
1294          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1295          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1296
1297          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1298          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1299
1300          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1301          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1302
1303          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1304          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1305
1306          IF ( i == 10 )  THEN
1307             EXIT
1308          ELSE
1309             i = i + 1
1310          ENDIF
1311
1312       ENDDO
1313
1314       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1315                          TRIM( gradients ), TRIM( slices )
1316    ENDIF
1317
1318
1319!
[1]1320!-- Listing of 1D-profiles
[151]1321    WRITE ( io, 325 )  dt_dopr_listing
[1353]1322    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1323       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1324    ENDIF
1325
1326!
1327!-- DATA output
1328    WRITE ( io, 330 )
[1353]1329    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1330       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1331    ENDIF
1332
1333!
1334!-- 1D-profiles
[346]1335    dopr_chr = 'Profile:'
[1]1336    IF ( dopr_n /= 0 )  THEN
1337       WRITE ( io, 331 )
1338
1339       output_format = ''
[1783]1340       output_format = netcdf_data_format_string
1341       IF ( netcdf_deflate == 0 )  THEN
1342          WRITE ( io, 344 )  output_format
1343       ELSE
1344          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1345       ENDIF
[1]1346
1347       DO  i = 1, dopr_n
1348          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
1349          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
1350             WRITE ( io, 332 )  dopr_chr
1351             dopr_chr = '       :'
1352          ENDIF
1353       ENDDO
1354
1355       IF ( dopr_chr /= '' )  THEN
1356          WRITE ( io, 332 )  dopr_chr
1357       ENDIF
1358       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
[1353]1359       IF ( skip_time_dopr /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dopr
[1]1360    ENDIF
1361
1362!
1363!-- 2D-arrays
1364    DO  av = 0, 1
1365
1366       i = 1
1367       do2d_xy = ''
1368       do2d_xz = ''
1369       do2d_yz = ''
1370       DO  WHILE ( do2d(av,i) /= ' ' )
1371
1372          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
1373          do2d_mode = do2d(av,i)(l-1:l)
1374
1375          SELECT CASE ( do2d_mode )
1376             CASE ( 'xy' )
1377                ll = LEN_TRIM( do2d_xy )
1378                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1379             CASE ( 'xz' )
1380                ll = LEN_TRIM( do2d_xz )
1381                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1382             CASE ( 'yz' )
1383                ll = LEN_TRIM( do2d_yz )
1384                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1385          END SELECT
1386
1387          i = i + 1
1388
1389       ENDDO
1390
1391       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
1392              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
[1327]1393              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) ) )  THEN
[1]1394
1395          IF (  av == 0 )  THEN
1396             WRITE ( io, 334 )  ''
1397          ELSE
1398             WRITE ( io, 334 )  '(time-averaged)'
1399          ENDIF
1400
1401          IF ( do2d_at_begin )  THEN
1402             begin_chr = 'and at the start'
1403          ELSE
1404             begin_chr = ''
1405          ENDIF
1406
1407          output_format = ''
[1783]1408          output_format = netcdf_data_format_string
1409          IF ( netcdf_deflate == 0 )  THEN
1410             WRITE ( io, 344 )  output_format
1411          ELSE
1412             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1413          ENDIF
[1]1414
1415          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
1416             i = 1
1417             slices = '/'
1418             coordinates = '/'
1419!
[1551]1420!--          Building strings with index and coordinate information of the
[1]1421!--          slices
1422             DO  WHILE ( section(i,1) /= -9999 )
1423
1424                WRITE (section_chr,'(I5)')  section(i,1)
1425                section_chr = ADJUSTL( section_chr )
1426                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1427
[206]1428                IF ( section(i,1) == -1 )  THEN
[1353]1429                   WRITE (coor_chr,'(F10.1)')  -1.0_wp
[206]1430                ELSE
1431                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
1432                ENDIF
[1]1433                coor_chr = ADJUSTL( coor_chr )
1434                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1435
1436                i = i + 1
1437             ENDDO
1438             IF ( av == 0 )  THEN
1439                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
1440                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
1441                                   TRIM( coordinates )
[1353]1442                IF ( skip_time_do2d_xy /= 0.0_wp )  THEN
[1]1443                   WRITE ( io, 339 )  skip_time_do2d_xy
1444                ENDIF
1445             ELSE
1446                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
1447                                   TRIM( begin_chr ), averaging_interval, &
1448                                   dt_averaging_input, 'k', TRIM( slices ), &
1449                                   TRIM( coordinates )
[1353]1450                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1451                   WRITE ( io, 339 )  skip_time_data_output_av
1452                ENDIF
1453             ENDIF
[1308]1454             IF ( netcdf_data_format > 4 )  THEN
1455                WRITE ( io, 352 )  ntdim_2d_xy(av)
1456             ELSE
1457                WRITE ( io, 353 )
1458             ENDIF
[1]1459          ENDIF
1460
1461          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
1462             i = 1
1463             slices = '/'
1464             coordinates = '/'
1465!
[1551]1466!--          Building strings with index and coordinate information of the
[1]1467!--          slices
1468             DO  WHILE ( section(i,2) /= -9999 )
1469
1470                WRITE (section_chr,'(I5)')  section(i,2)
1471                section_chr = ADJUSTL( section_chr )
1472                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1473
1474                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
1475                coor_chr = ADJUSTL( coor_chr )
1476                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1477
1478                i = i + 1
1479             ENDDO
1480             IF ( av == 0 )  THEN
1481                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
1482                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
1483                                   TRIM( coordinates )
[1353]1484                IF ( skip_time_do2d_xz /= 0.0_wp )  THEN
[1]1485                   WRITE ( io, 339 )  skip_time_do2d_xz
1486                ENDIF
1487             ELSE
1488                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
1489                                   TRIM( begin_chr ), averaging_interval, &
1490                                   dt_averaging_input, 'j', TRIM( slices ), &
1491                                   TRIM( coordinates )
[1353]1492                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1493                   WRITE ( io, 339 )  skip_time_data_output_av
1494                ENDIF
1495             ENDIF
[1308]1496             IF ( netcdf_data_format > 4 )  THEN
1497                WRITE ( io, 352 )  ntdim_2d_xz(av)
1498             ELSE
1499                WRITE ( io, 353 )
1500             ENDIF
[1]1501          ENDIF
1502
1503          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
1504             i = 1
1505             slices = '/'
1506             coordinates = '/'
1507!
[1551]1508!--          Building strings with index and coordinate information of the
[1]1509!--          slices
1510             DO  WHILE ( section(i,3) /= -9999 )
1511
1512                WRITE (section_chr,'(I5)')  section(i,3)
1513                section_chr = ADJUSTL( section_chr )
1514                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1515
1516                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
1517                coor_chr = ADJUSTL( coor_chr )
1518                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1519
1520                i = i + 1
1521             ENDDO
1522             IF ( av == 0 )  THEN
1523                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
1524                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
1525                                   TRIM( coordinates )
[1353]1526                IF ( skip_time_do2d_yz /= 0.0_wp )  THEN
[1]1527                   WRITE ( io, 339 )  skip_time_do2d_yz
1528                ENDIF
1529             ELSE
1530                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
1531                                   TRIM( begin_chr ), averaging_interval, &
1532                                   dt_averaging_input, 'i', TRIM( slices ), &
1533                                   TRIM( coordinates )
[1353]1534                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1535                   WRITE ( io, 339 )  skip_time_data_output_av
1536                ENDIF
1537             ENDIF
[1308]1538             IF ( netcdf_data_format > 4 )  THEN
1539                WRITE ( io, 352 )  ntdim_2d_yz(av)
1540             ELSE
1541                WRITE ( io, 353 )
1542             ENDIF
[1]1543          ENDIF
1544
1545       ENDIF
1546
1547    ENDDO
1548
1549!
1550!-- 3d-arrays
1551    DO  av = 0, 1
1552
1553       i = 1
1554       do3d_chr = ''
1555       DO  WHILE ( do3d(av,i) /= ' ' )
1556
1557          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
1558          i = i + 1
1559
1560       ENDDO
1561
1562       IF ( do3d_chr /= '' )  THEN
1563          IF ( av == 0 )  THEN
1564             WRITE ( io, 336 )  ''
1565          ELSE
1566             WRITE ( io, 336 )  '(time-averaged)'
1567          ENDIF
1568
[1783]1569          output_format = netcdf_data_format_string
1570          IF ( netcdf_deflate == 0 )  THEN
1571             WRITE ( io, 344 )  output_format
1572          ELSE
1573             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1574          ENDIF
[1]1575
1576          IF ( do3d_at_begin )  THEN
1577             begin_chr = 'and at the start'
1578          ELSE
1579             begin_chr = ''
1580          ENDIF
1581          IF ( av == 0 )  THEN
1582             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
1583                                zu(nz_do3d), nz_do3d
1584          ELSE
1585             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
1586                                TRIM( begin_chr ), averaging_interval, &
1587                                dt_averaging_input, zu(nz_do3d), nz_do3d
1588          ENDIF
1589
[1308]1590          IF ( netcdf_data_format > 4 )  THEN
1591             WRITE ( io, 352 )  ntdim_3d(av)
1592          ELSE
1593             WRITE ( io, 353 )
1594          ENDIF
1595
[1]1596          IF ( av == 0 )  THEN
[1353]1597             IF ( skip_time_do3d /= 0.0_wp )  THEN
[1]1598                WRITE ( io, 339 )  skip_time_do3d
1599             ENDIF
1600          ELSE
[1353]1601             IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1602                WRITE ( io, 339 )  skip_time_data_output_av
1603             ENDIF
1604          ENDIF
1605
1606       ENDIF
1607
1608    ENDDO
1609
1610!
[410]1611!-- masked arrays
1612    IF ( masks > 0 )  WRITE ( io, 345 )  &
1613         mask_scale_x, mask_scale_y, mask_scale_z
1614    DO  mid = 1, masks
1615       DO  av = 0, 1
1616
1617          i = 1
1618          domask_chr = ''
1619          DO  WHILE ( domask(mid,av,i) /= ' ' )
1620             domask_chr = TRIM( domask_chr ) // ' ' //  &
1621                          TRIM( domask(mid,av,i) ) // ','
1622             i = i + 1
1623          ENDDO
1624
1625          IF ( domask_chr /= '' )  THEN
1626             IF ( av == 0 )  THEN
1627                WRITE ( io, 346 )  '', mid
1628             ELSE
1629                WRITE ( io, 346 )  ' (time-averaged)', mid
1630             ENDIF
1631
[1783]1632             output_format = netcdf_data_format_string
[1308]1633!--          Parallel output not implemented for mask data, hence
1634!--          output_format must be adjusted.
1635             IF ( netcdf_data_format == 5 ) output_format = 'netCDF4/HDF5'
1636             IF ( netcdf_data_format == 6 ) output_format = 'netCDF4/HDF5 classic'
[1783]1637             IF ( netcdf_deflate == 0 )  THEN
1638                WRITE ( io, 344 )  output_format
1639             ELSE
1640                WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1641             ENDIF
[410]1642
1643             IF ( av == 0 )  THEN
1644                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1645             ELSE
1646                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1647                                   averaging_interval, dt_averaging_input
1648             ENDIF
1649
1650             IF ( av == 0 )  THEN
[1353]1651                IF ( skip_time_domask(mid) /= 0.0_wp )  THEN
[410]1652                   WRITE ( io, 339 )  skip_time_domask(mid)
1653                ENDIF
1654             ELSE
[1353]1655                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[410]1656                   WRITE ( io, 339 )  skip_time_data_output_av
1657                ENDIF
1658             ENDIF
1659!
1660!--          output locations
1661             DO  dim = 1, 3
[1353]1662                IF ( mask(mid,dim,1) >= 0.0_wp )  THEN
[410]1663                   count = 0
[1353]1664                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0_wp )
[410]1665                      count = count + 1
1666                   ENDDO
1667                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1668                                      mask(mid,dim,:count)
[1353]1669                ELSEIF ( mask_loop(mid,dim,1) < 0.0_wp .AND.  &
1670                         mask_loop(mid,dim,2) < 0.0_wp .AND.  &
1671                         mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1672                   WRITE ( io, 350 )  dir(dim), dir(dim)
[1353]1673                ELSEIF ( mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1674                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1675                                      mask_loop(mid,dim,1:2)
1676                ELSE
1677                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1678                                      mask_loop(mid,dim,1:3)
1679                ENDIF
1680             ENDDO
1681          ENDIF
1682
1683       ENDDO
1684    ENDDO
1685
1686!
[1]1687!-- Timeseries
[1322]1688    IF ( dt_dots /= 9999999.9_wp )  THEN
[1]1689       WRITE ( io, 340 )
1690
[1783]1691       output_format = netcdf_data_format_string
1692       IF ( netcdf_deflate == 0 )  THEN
1693          WRITE ( io, 344 )  output_format
1694       ELSE
1695          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1696       ENDIF
[1]1697       WRITE ( io, 341 )  dt_dots
1698    ENDIF
1699
1700#if defined( __dvrp_graphics )
1701!
1702!-- Dvrp-output
[1322]1703    IF ( dt_dvrp /= 9999999.9_wp )  THEN
[1]1704       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1705                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1706       i = 1
1707       l = 0
[336]1708       m = 0
[1]1709       DO WHILE ( mode_dvrp(i) /= ' ' )
1710          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
[130]1711             READ ( mode_dvrp(i), '(10X,I2)' )  j
[1]1712             l = l + 1
1713             IF ( do3d(0,j) /= ' ' )  THEN
[336]1714                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1715                                   isosurface_color(:,l)
[1]1716             ENDIF
1717          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
[130]1718             READ ( mode_dvrp(i), '(6X,I2)' )  j
[336]1719             m = m + 1
1720             IF ( do2d(0,j) /= ' ' )  THEN
1721                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1722                                   slicer_range_limits_dvrp(:,m)
1723             ENDIF
[1]1724          ENDIF
1725          i = i + 1
1726       ENDDO
[237]1727
[336]1728       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1729                          superelevation_y, superelevation, clip_dvrp_l, &
1730                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1731
1732       IF ( TRIM( topography ) /= 'flat' )  THEN
1733          WRITE ( io, 366 )  topography_color
1734          IF ( cluster_size > 1 )  THEN
1735             WRITE ( io, 367 )  cluster_size
1736          ENDIF
[237]1737       ENDIF
1738
[1]1739    ENDIF
1740#endif
[1957]1741!
1742!-- Output of virtual flight information
1743    IF ( virtual_flight )  CALL flight_header( io )
[1]1744
1745!
[1833]1746!-- Output of spectra related quantities
1747    IF ( calculate_spectra )  CALL spectra_header( io )
[1]1748
1749    WRITE ( io, 99 )
1750
1751!
1752!-- Physical quantities
1753    WRITE ( io, 400 )
1754
1755!
1756!-- Geostrophic parameters
[1551]1757    WRITE ( io, 410 )  phi, omega, f, fs
[1]1758
1759!
1760!-- Other quantities
1761    WRITE ( io, 411 )  g
[1551]1762
[1179]1763    WRITE ( io, 412 )  TRIM( reference_state )
1764    IF ( use_single_reference_value )  THEN
[97]1765       IF ( ocean )  THEN
[1179]1766          WRITE ( io, 413 )  prho_reference
[97]1767       ELSE
[1179]1768          WRITE ( io, 414 )  pt_reference
[97]1769       ENDIF
1770    ENDIF
[1]1771
1772!
1773!-- Cloud physics parameters
[1299]1774    IF ( cloud_physics )  THEN
[57]1775       WRITE ( io, 415 )
1776       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
[1822]1777       IF ( microphysics_seifert )  THEN
[1353]1778          WRITE ( io, 510 ) 1.0E-6_wp * nc_const
[1822]1779          WRITE ( io, 511 ) c_sedimentation
[1115]1780       ENDIF
[1]1781    ENDIF
1782
1783!
[824]1784!-- Cloud physcis parameters / quantities / numerical methods
1785    WRITE ( io, 430 )
1786    IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets)  THEN
1787       WRITE ( io, 431 )
1788    ELSEIF ( humidity  .AND.  cloud_physics )  THEN
1789       WRITE ( io, 432 )
[1496]1790       IF ( cloud_top_radiation )  WRITE ( io, 132 )
[1822]1791       IF ( microphysics_kessler )  THEN
1792          WRITE ( io, 133 )
1793       ELSEIF ( microphysics_seifert )  THEN
[1831]1794          IF ( cloud_water_sedimentation )  WRITE ( io, 506 )
[1822]1795          WRITE ( io, 505 )
[1831]1796          IF ( collision_turbulence )  WRITE ( io, 507 )
[1822]1797          IF ( ventilation_effect )  WRITE ( io, 508 )
1798          IF ( limiter_sedimentation )  WRITE ( io, 509 )
[1115]1799       ENDIF
[824]1800    ELSEIF ( humidity  .AND.  cloud_droplets )  THEN
1801       WRITE ( io, 433 )
1802       IF ( curvature_solution_effects )  WRITE ( io, 434 )
[825]1803       IF ( collision_kernel /= 'none' )  THEN
1804          WRITE ( io, 435 )  TRIM( collision_kernel )
[828]1805          IF ( collision_kernel(6:9) == 'fast' )  THEN
1806             WRITE ( io, 436 )  radius_classes, dissipation_classes
1807          ENDIF
[825]1808       ELSE
[828]1809          WRITE ( io, 437 )
[825]1810       ENDIF
[824]1811    ENDIF
1812
1813!
[1]1814!-- LES / turbulence parameters
1815    WRITE ( io, 450 )
1816
1817!--
1818! ... LES-constants used must still be added here
1819!--
1820    IF ( constant_diffusion )  THEN
1821       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1822                          prandtl_number
1823    ENDIF
1824    IF ( .NOT. constant_diffusion)  THEN
[1353]1825       IF ( e_init > 0.0_wp )  WRITE ( io, 455 )  e_init
1826       IF ( e_min > 0.0_wp )  WRITE ( io, 454 )  e_min
[1]1827       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1828    ENDIF
1829
1830!
1831!-- Special actions during the run
1832    WRITE ( io, 470 )
1833    IF ( create_disturbances )  THEN
1834       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1835                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1836                          zu(disturbance_level_ind_t), disturbance_level_ind_t
[707]1837       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1]1838          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1839       ELSE
1840          WRITE ( io, 473 )  disturbance_energy_limit
1841       ENDIF
1842       WRITE ( io, 474 )  TRIM( random_generator )
1843    ENDIF
[1353]1844    IF ( pt_surface_initial_change /= 0.0_wp )  THEN
[1]1845       WRITE ( io, 475 )  pt_surface_initial_change
1846    ENDIF
[1353]1847    IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1848       WRITE ( io, 476 )  q_surface_initial_change       
1849    ENDIF
[1353]1850    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1851       WRITE ( io, 477 )  q_surface_initial_change       
1852    ENDIF
1853
[60]1854    IF ( particle_advection )  THEN
[1]1855!
[60]1856!--    Particle attributes
1857       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1858                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
[1359]1859                          end_time_prel
[60]1860       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1861       IF ( random_start_position )  WRITE ( io, 481 )
[1575]1862       IF ( seed_follows_topography )  WRITE ( io, 496 )
[60]1863       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1864       WRITE ( io, 495 )  total_number_of_particles
[1322]1865       IF ( dt_write_particle_data /= 9999999.9_wp )  THEN
[60]1866          WRITE ( io, 485 )  dt_write_particle_data
[1327]1867          IF ( netcdf_data_format > 1 )  THEN
1868             output_format = 'netcdf (64 bit offset) and binary'
[1]1869          ELSE
[1327]1870             output_format = 'netcdf and binary'
[1]1871          ENDIF
[1783]1872          IF ( netcdf_deflate == 0 )  THEN
1873             WRITE ( io, 344 )  output_format
1874          ELSE
1875             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1876          ENDIF
[1]1877       ENDIF
[1322]1878       IF ( dt_dopts /= 9999999.9_wp )  WRITE ( io, 494 )  dt_dopts
[60]1879       IF ( write_particle_statistics )  WRITE ( io, 486 )
[1]1880
[60]1881       WRITE ( io, 487 )  number_of_particle_groups
[1]1882
[60]1883       DO  i = 1, number_of_particle_groups
[1322]1884          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9_wp )  THEN
[1353]1885             WRITE ( io, 490 )  i, 0.0_wp
[60]1886             WRITE ( io, 492 )
[1]1887          ELSE
[60]1888             WRITE ( io, 490 )  i, radius(i)
[1353]1889             IF ( density_ratio(i) /= 0.0_wp )  THEN
[60]1890                WRITE ( io, 491 )  density_ratio(i)
1891             ELSE
1892                WRITE ( io, 492 )
1893             ENDIF
[1]1894          ENDIF
[60]1895          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1896                             pdx(i), pdy(i), pdz(i)
[336]1897          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
[60]1898       ENDDO
[1]1899
[60]1900    ENDIF
[1]1901
[60]1902
[1]1903!
1904!-- Parameters of 1D-model
1905    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1906       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1907                          mixing_length_1d, dissipation_1d
1908       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1909          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1910       ENDIF
1911    ENDIF
1912
1913!
[1551]1914!-- User-defined information
[1]1915    CALL user_header( io )
1916
1917    WRITE ( io, 99 )
1918
1919!
1920!-- Write buffer contents to disc immediately
[1808]1921    FLUSH( io )
[1]1922
1923!
1924!-- Here the FORMATs start
1925
1926 99 FORMAT (1X,78('-'))
[1468]1927100 FORMAT (/1X,'******************************',4X,44('-')/        &
1928            1X,'* ',A,' *',4X,A/                               &
1929            1X,'******************************',4X,44('-'))
1930101 FORMAT (35X,'coupled run using MPI-',I1,': ',A/ &
1931            35X,42('-'))
1932102 FORMAT (/' Date:                 ',A8,4X,'Run:       ',A20/      &
1933            ' Time:                 ',A8,4X,'Run-No.:   ',I2.2/     &
[1106]1934            ' Run on host:        ',A10)
[1]1935#if defined( __parallel )
[1468]1936103 FORMAT (' Number of PEs:',10X,I6,4X,'Processor grid (x,y): (',I4,',',I4, &
[1]1937              ')',1X,A)
[1468]1938104 FORMAT (' Number of PEs:',10X,I6,4X,'Tasks:',I4,'   threads per task:',I4/ &
1939              35X,'Processor grid (x,y): (',I4,',',I4,')',1X,A)
1940105 FORMAT (35X,'One additional PE is used to handle'/37X,'the dvrp output!')
1941106 FORMAT (35X,'A 1d-decomposition along x is forced'/ &
1942            35X,'because the job is running on an SMP-cluster')
1943107 FORMAT (35X,'A 1d-decomposition along ',A,' is used')
1944108 FORMAT (35X,'Max. # of parallel I/O streams is ',I5)
1945109 FORMAT (35X,'Precursor run for coupled atmos-ocean run'/ &
1946            35X,42('-'))
1947114 FORMAT (35X,'Coupled atmosphere-ocean run following'/ &
1948            35X,'independent precursor runs'/             &
1949            35X,42('-'))
[1]1950#endif
1951110 FORMAT (/' Numerical Schemes:'/ &
1952             ' -----------------'/)
[2037]1953121 FORMAT (' --> Use the ',A,' approximation for the model equations.')
[1]1954111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1955112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
[1697]1956            '     Iterations (initial/other): ',I3,'/',I3,'  omega =',F6.3)
[1]1957113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1958                  ' or Upstream')
[1216]1959115 FORMAT ('     FFT and transpositions are overlapping')
[1]1960116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1961                  ' or Upstream')
1962118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
[1106]1963119 FORMAT (' --> Galilei-Transform applied to horizontal advection:'/ &
1964            '     translation velocity = ',A/ &
[1]1965            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1966122 FORMAT (' --> Time differencing scheme: ',A)
[108]1967123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
[1697]1968            '     maximum damping coefficient:',F6.3, ' 1/s')
[1]1969129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1970130 FORMAT (' --> Additional prognostic equation for the total water content')
[940]1971131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
1972                  F6.2, ' K assumed')
[824]1973132 FORMAT ('     Parameterization of long-wave radiation processes via'/ &
[1]1974            '     effective emissivity scheme')
[824]1975133 FORMAT ('     Precipitation parameterization via Kessler-Scheme')
[1]1976134 FORMAT (' --> Additional prognostic equation for a passive scalar')
[1575]1977135 FORMAT (' --> Solve perturbation pressure via ',A,' method (', &
[1]1978                  A,'-cycle)'/ &
1979            '     number of grid levels:                   ',I2/ &
1980            '     Gauss-Seidel red/black iterations:       ',I2)
1981136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1982                  I3,')')
1983137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1984            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1985                  I3,')'/ &
1986            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1987                  I3,')')
[63]1988139 FORMAT (' --> Loop optimization method: ',A)
[1]1989140 FORMAT ('     maximum residual allowed:                ',E10.3)
1990141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1991142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1992                  'step')
[87]1993143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1994                  'kinetic energy')
[927]1995144 FORMAT ('     masking method is used')
[1]1996150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
[241]1997                  'conserved'/ &
1998            '     using the ',A,' mode')
1999151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
[306]2000152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
2001           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
2002           /'     starting from dp_level_b =', F8.3, 'm', A /)
[1365]2003160 FORMAT (//' Large scale forcing and nudging:'/ &
2004              ' -------------------------------'/)
2005161 FORMAT (' --> No large scale forcing from external is used (default) ')
2006162 FORMAT (' --> Large scale forcing from external file LSF_DATA is used: ')
2007163 FORMAT ('     - large scale advection tendencies ')
2008164 FORMAT ('     - large scale subsidence velocity w_subs ')
2009165 FORMAT ('     - large scale subsidence tendencies ')
2010167 FORMAT ('     - and geostrophic wind components ug and vg')
2011168 FORMAT (' --> Large-scale vertical motion is used in the ', &
[1299]2012                  'prognostic equation(s) for')
[1365]2013169 FORMAT ('     the scalar(s) only')
2014170 FORMAT (' --> Nudging is used')
2015171 FORMAT (' --> No nudging is used (default) ')
2016180 FORMAT ('     - prescribed surface values for temperature')
[1376]2017181 FORMAT ('     - prescribed surface fluxes for temperature')
2018182 FORMAT ('     - prescribed surface values for humidity')
[1365]2019183 FORMAT ('     - prescribed surface fluxes for humidity')
[1]2020200 FORMAT (//' Run time and time step information:'/ &
2021             ' ----------------------------------'/)
[1106]2022201 FORMAT ( ' Timestep:             variable     maximum value: ',F6.3,' s', &
[1697]2023             '    CFL-factor:',F5.2)
[1106]2024202 FORMAT ( ' Timestep:          dt = ',F6.3,' s'/)
2025203 FORMAT ( ' Start time:          ',F9.3,' s'/ &
2026             ' End time:            ',F9.3,' s')
[1]2027204 FORMAT ( A,F9.3,' s')
2028205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
[1106]2029206 FORMAT (/' Time reached:        ',F9.3,' s'/ &
2030             ' CPU-time used:       ',F9.3,' s     per timestep:               ', &
2031               '  ',F9.3,' s'/                                                    &
[1111]2032             '                                      per second of simulated tim', &
[1]2033               'e: ',F9.3,' s')
[1106]2034207 FORMAT ( ' Coupling start time: ',F9.3,' s')
[1]2035250 FORMAT (//' Computational grid and domain size:'/ &
2036              ' ----------------------------------'// &
2037              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
2038              ' m    dz =    ',F7.3,' m'/ &
2039              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
2040              ' m  z(u) = ',F10.3,' m'/)
2041252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
[1697]2042              ' factor:',F6.3/ &
[1]2043            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
2044254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
2045            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
2046260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
2047             ' degrees')
[1551]2048270 FORMAT (//' Topography information:'/ &
2049              ' ----------------------'// &
[1]2050              1X,'Topography: ',A)
2051271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
2052              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
2053                ' / ',I4)
[240]2054272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
2055              ' direction' / &
2056              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
2057              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
[2232]2058273 FORMAT (  ' Tunnel of infinite length in ',A, &
2059              ' direction' / &
2060              ' Tunnel height: ', F6.2, / &
2061              ' Tunnel-wall depth: ', F6.2      / &
2062              ' Tunnel width: ', F6.2 )
2063274 FORMAT (  ' Tunnel in ', A, ' direction.' / &
2064              ' Tunnel height: ', F6.2, / &   
2065              ' Tunnel-wall depth: ', F6.2      / &
2066              ' Tunnel width: ', F6.2, / &
2067              ' Tunnel length: ', F6.2 )
[256]2068278 FORMAT (' Topography grid definition convention:'/ &
2069            ' cell edge (staggered grid points'/  &
2070            ' (u in x-direction, v in y-direction))' /)
2071279 FORMAT (' Topography grid definition convention:'/ &
2072            ' cell center (scalar grid points)' /)
[1]2073300 FORMAT (//' Boundary conditions:'/ &
2074             ' -------------------'// &
2075             '                     p                    uv             ', &
[1551]2076             '                     pt'// &
[1]2077             ' B. bound.: ',A/ &
2078             ' T. bound.: ',A)
[97]2079301 FORMAT (/'                     ',A// &
[1]2080             ' B. bound.: ',A/ &
2081             ' T. bound.: ',A)
[19]2082303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
2083304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
2084305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
2085               'computational u,v-level:'// &
[1697]2086             '       zp = ',F6.2,' m   z0 =',F7.4,' m   z0h =',F8.5,&
2087             ' m   kappa =',F5.2/ &
2088             '       Rif value range:   ',F8.2,' <= rif <=',F6.2)
[97]2089306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
[1]2090307 FORMAT ('       Heatflux has a random normal distribution')
2091308 FORMAT ('       Predefined surface temperature')
[97]2092309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
[1]2093310 FORMAT (//'    1D-Model:'// &
2094             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
[1960]2095311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' kg/kg m/s')
[1]2096312 FORMAT ('       Predefined surface humidity')
2097313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
2098314 FORMAT ('       Predefined scalar value at the surface')
[1992]2099302 FORMAT ('       Predefined constant scalarflux:   ',F9.6,' kg/(m**2 s)')
2100315 FORMAT ('       Humidity flux at top surface is 0.0')
[102]2101316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
2102                    'atmosphere model')
[1]2103317 FORMAT (//' Lateral boundaries:'/ &
2104            '       left/right:  ',A/    &
2105            '       north/south: ',A)
[1159]2106318 FORMAT (/'       use_cmax: ',L1 / &
2107            '       pt damping layer width = ',F8.2,' m, pt ', &
[1697]2108                    'damping factor =',F7.4)
[151]2109319 FORMAT ('       turbulence recycling at inflow switched on'/ &
2110            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
2111            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
2112320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
[103]2113            '                                          v: ',F9.6,' m**2/s**2')
[1365]2114321 FORMAT (//' Initial profiles:'/ &
2115              ' ----------------')
[1560]2116322 FORMAT ('       turbulence recycling at inflow switched on'/ &
2117            '       y shift of the recycled inflow turbulence switched on'/ &
2118            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
[1592]2119            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m'/)
[2050]2120323 FORMAT ('       turbulent outflow conditon switched on'/ &
2121            '       position of outflow source plane: ',F7.1,' m   ', &
2122                    'grid index: ', I4)
[151]2123325 FORMAT (//' List output:'/ &
[1]2124             ' -----------'//  &
2125            '    1D-Profiles:'/    &
2126            '       Output every             ',F8.2,' s')
[151]2127326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
[1]2128            '       Averaging input every    ',F8.2,' s')
2129330 FORMAT (//' Data output:'/ &
2130             ' -----------'/)
2131331 FORMAT (/'    1D-Profiles:')
2132332 FORMAT (/'       ',A)
2133333 FORMAT ('       Output every             ',F8.2,' s',/ &
2134            '       Time averaged over       ',F8.2,' s'/ &
2135            '       Averaging input every    ',F8.2,' s')
2136334 FORMAT (/'    2D-Arrays',A,':')
2137335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2138            '       Output every             ',F8.2,' s  ',A/ &
2139            '       Cross sections at ',A1,' = ',A/ &
2140            '       scalar-coordinates:   ',A,' m'/)
2141336 FORMAT (/'    3D-Arrays',A,':')
2142337 FORMAT (/'       Arrays: ',A/ &
2143            '       Output every             ',F8.2,' s  ',A/ &
2144            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
2145339 FORMAT ('       No output during initial ',F8.2,' s')
2146340 FORMAT (/'    Time series:')
2147341 FORMAT ('       Output every             ',F8.2,' s'/)
2148342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2149            '       Output every             ',F8.2,' s  ',A/ &
2150            '       Time averaged over       ',F8.2,' s'/ &
2151            '       Averaging input every    ',F8.2,' s'/ &
2152            '       Cross sections at ',A1,' = ',A/ &
2153            '       scalar-coordinates:   ',A,' m'/)
2154343 FORMAT (/'       Arrays: ',A/ &
2155            '       Output every             ',F8.2,' s  ',A/ &
2156            '       Time averaged over       ',F8.2,' s'/ &
2157            '       Averaging input every    ',F8.2,' s'/ &
2158            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
[292]2159344 FORMAT ('       Output format: ',A/)
[410]2160345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
2161            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
2162            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
2163            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
2164346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
2165347 FORMAT ('       Variables: ',A/ &
2166            '       Output every             ',F8.2,' s')
2167348 FORMAT ('       Variables: ',A/ &
2168            '       Output every             ',F8.2,' s'/ &
2169            '       Time averaged over       ',F8.2,' s'/ &
2170            '       Averaging input every    ',F8.2,' s')
2171349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2172            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
2173            13('       ',8(F8.2,',')/) )
2174350 FORMAT (/'       Output locations in ',A,'-direction: ', &
2175            'all gridpoints along ',A,'-direction (default).' )
2176351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2177            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
2178            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
[1313]2179352 FORMAT  (/'       Number of output time levels allowed: ',I3 /)
2180353 FORMAT  (/'       Number of output time levels allowed: unlimited' /)
[1783]2181354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
[1]2182#if defined( __dvrp_graphics )
2183360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
2184            '       Output every      ',F7.1,' s'/ &
2185            '       Output mode:      ',A/ &
2186            '       Host / User:      ',A,' / ',A/ &
2187            '       Directory:        ',A// &
2188            '       The sequence contains:')
[337]2189361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
2190            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
2191362 FORMAT (/'       Slicer plane ',A/ &
[336]2192            '       Slicer limits: [',F6.2,',',F6.2,']')
[337]2193365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
[336]2194            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
2195                     ')'/ &
2196            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
2197            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
[337]2198366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
[336]2199367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
[1]2200#endif
2201400 FORMAT (//' Physical quantities:'/ &
2202              ' -------------------'/)
[1551]2203410 FORMAT ('    Geograph. latitude  :   phi    = ',F4.1,' degr'/   &
[1697]2204            '    Angular velocity    :   omega  =',E10.3,' rad/s'/  &
[1551]2205            '    Coriolis parameter  :   f      = ',F9.6,' 1/s'/    &
2206            '                            f*     = ',F9.6,' 1/s')
2207411 FORMAT (/'    Gravity             :   g      = ',F4.1,' m/s**2')
[1179]2208412 FORMAT (/'    Reference state used in buoyancy terms: ',A)
2209413 FORMAT ('       Reference density in buoyancy terms: ',F8.3,' kg/m**3')
2210414 FORMAT ('       Reference temperature in buoyancy terms: ',F8.4,' K')
[1551]2211415 FORMAT (/' Cloud physics parameters:'/ &
2212             ' ------------------------'/)
2213416 FORMAT ('    Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
2214            '    Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
[1697]2215            '    Density of air     :   rho_0 =',F6.3,' kg/m**3'/  &
[1551]2216            '    Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
[1697]2217            '    Vapourization heat :   L_v   =',E9.2,' J/kg')
[1551]2218417 FORMAT ('    Geograph. longitude :   lambda = ',F4.1,' degr')
2219418 FORMAT (/'    Day of the year at model start :   day_init      =     ',I3 &
2220            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
[1]2221420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
2222            '       Height:        ',A,'  m'/ &
2223            '       Temperature:   ',A,'  K'/ &
2224            '       Gradient:      ',A,'  K/100m'/ &
2225            '       Gridpoint:     ',A)
2226421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
2227            '       Height:      ',A,'  m'/ &
2228            '       Humidity:    ',A,'  kg/kg'/ &
2229            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
2230            '       Gridpoint:   ',A)
2231422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
2232            '       Height:                  ',A,'  m'/ &
2233            '       Scalar concentration:    ',A,'  kg/m**3'/ &
2234            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
2235            '       Gridpoint:               ',A)
2236423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
2237            '       Height:      ',A,'  m'/ &
2238            '       ug:          ',A,'  m/s'/ &
2239            '       Gradient:    ',A,'  1/100s'/ &
2240            '       Gridpoint:   ',A)
2241424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
2242            '       Height:      ',A,'  m'/ &
[97]2243            '       vg:          ',A,'  m/s'/ &
[1]2244            '       Gradient:    ',A,'  1/100s'/ &
2245            '       Gridpoint:   ',A)
[97]2246425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
2247            '       Height:     ',A,'  m'/ &
2248            '       Salinity:   ',A,'  psu'/ &
2249            '       Gradient:   ',A,'  psu/100m'/ &
2250            '       Gridpoint:  ',A)
[411]2251426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
2252            '       Height:      ',A,'  m'/ &
2253            '       w_subs:      ',A,'  m/s'/ &
2254            '       Gradient:    ',A,'  (m/s)/100m'/ &
2255            '       Gridpoint:   ',A)
[767]2256427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
2257                  ' profiles')
[1241]2258428 FORMAT (/'    Initial profiles (u, v, pt, q) are taken from file '/ &
2259             '    NUDGING_DATA')
[824]2260430 FORMAT (//' Cloud physics quantities / methods:'/ &
2261              ' ----------------------------------'/)
[1960]2262431 FORMAT ('    Humidity is considered, bu no condensation')
[824]2263432 FORMAT ('    Bulk scheme with liquid water potential temperature and'/ &
2264            '    total water content is used.'/ &
2265            '    Condensation is parameterized via 0% - or 100% scheme.')
2266433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
2267                 'icle model')
2268434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
2269                 ' droplets < 1.0E-6 m')
[825]2270435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
[828]2271436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
2272                    'are used'/ &
2273            '          number of radius classes:       ',I3,'    interval ', &
2274                       '[1.0E-6,2.0E-4] m'/ &
2275            '          number of dissipation classes:   ',I2,'    interval ', &
2276                       '[0,1000] cm**2/s**3')
2277437 FORMAT ('    Droplet collision is switched off')
[1]2278450 FORMAT (//' LES / Turbulence quantities:'/ &
2279              ' ---------------------------'/)
[824]2280451 FORMAT ('    Diffusion coefficients are constant:'/ &
2281            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
[1697]2282453 FORMAT ('    Mixing length is limited to',F5.2,' * z')
[824]2283454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
2284455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
[1]2285470 FORMAT (//' Actions during the simulation:'/ &
2286              ' -----------------------------'/)
[94]2287471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
[1697]2288            '    Disturbance amplitude           :    ',F5.2, ' m/s'/       &
[94]2289            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
2290            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
[1]2291472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
2292                 ' to i/j =',I4)
2293473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
[1697]2294                 F6.3, ' m**2/s**2')
[1]2295474 FORMAT ('    Random number generator used    : ',A/)
2296475 FORMAT ('    The surface temperature is increased (or decreased, ', &
2297                 'respectively, if'/ &
2298            '    the value is negative) by ',F5.2,' K at the beginning of the',&
2299                 ' 3D-simulation'/)
2300476 FORMAT ('    The surface humidity is increased (or decreased, ',&
2301                 'respectively, if the'/ &
2302            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
2303                 ' the 3D-simulation'/)
2304477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
2305                 'respectively, if the'/ &
2306            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
2307                 ' the 3D-simulation'/)
2308480 FORMAT ('    Particles:'/ &
2309            '    ---------'// &
2310            '       Particle advection is active (switched on at t = ', F7.1, &
2311                    ' s)'/ &
2312            '       Start of new particle generations every  ',F6.1,' s'/ &
2313            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
2314            '                            bottom:     ', A, ' top:         ', A/&
2315            '       Maximum particle age:                 ',F9.1,' s'/ &
[1359]2316            '       Advection stopped at t = ',F9.1,' s'/)
[1]2317481 FORMAT ('       Particles have random start positions'/)
[336]2318482 FORMAT ('          Particles are advected only horizontally'/)
[1]2319485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
2320486 FORMAT ('       Particle statistics are written on file'/)
2321487 FORMAT ('       Number of particle groups: ',I2/)
2322488 FORMAT ('       SGS velocity components are used for particle advection'/ &
[1697]2323            '          minimum timestep for advection:', F8.5/)
[1]2324489 FORMAT ('       Number of particles simultaneously released at each ', &
2325                    'point: ', I5/)
2326490 FORMAT ('       Particle group ',I2,':'/ &
2327            '          Particle radius: ',E10.3, 'm')
2328491 FORMAT ('          Particle inertia is activated'/ &
[1697]2329            '             density_ratio (rho_fluid/rho_particle) =',F6.3/)
[1]2330492 FORMAT ('          Particles are advected only passively (no inertia)'/)
2331493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
2332            '                                         y:',F8.1,' - ',F8.1,' m'/&
2333            '                                         z:',F8.1,' - ',F8.1,' m'/&
2334            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
2335                       ' m  dz = ',F8.1,' m'/)
2336494 FORMAT ('       Output of particle time series in NetCDF format every ', &
2337                    F8.2,' s'/)
2338495 FORMAT ('       Number of particles in total domain: ',I10/)
[1575]2339496 FORMAT ('       Initial vertical particle positions are interpreted ', &
2340                    'as relative to the given topography')
[1]2341500 FORMAT (//' 1D-Model parameters:'/                           &
2342              ' -------------------'//                           &
2343            '    Simulation time:                   ',F8.1,' s'/ &
2344            '    Run-controll output every:         ',F8.1,' s'/ &
2345            '    Vertical profile output every:     ',F8.1,' s'/ &
2346            '    Mixing length calculation:         ',A/         &
2347            '    Dissipation calculation:           ',A/)
2348502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
[667]2349503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2350504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
[1115]2351505 FORMAT ('    Precipitation parameterization via Seifert-Beheng-Scheme')
[1831]2352506 FORMAT ('    Cloud water sedimentation parameterization via Stokes law')
[1115]2353507 FORMAT ('    Turbulence effects on precipitation process')
2354508 FORMAT ('    Ventilation effects on evaporation of rain drops')
2355509 FORMAT ('    Slope limiter used for sedimentation process')
[1551]2356510 FORMAT ('    Droplet density    :   N_c   = ',F6.1,' 1/cm**3')
2357511 FORMAT ('    Sedimentation Courant number:                  '/&
[1697]2358            '                               C_s   =',F4.1,'        ')
[1429]2359512 FORMAT (/' Date:                 ',A8,6X,'Run:       ',A20/      &
2360            ' Time:                 ',A8,6X,'Run-No.:   ',I2.2/     &
2361            ' Run on host:        ',A10,6X,'En-No.:    ',I2.2)
[1791]2362600 FORMAT (/' Nesting informations:'/ &
2363            ' --------------------'/ &
[1797]2364            ' Nesting mode:                     ',A/ &
2365            ' Nesting-datatransfer mode:        ',A// &
[1791]2366            ' Nest id  parent  number   lower left coordinates   name'/ &
2367            ' (*=me)     id    of PEs      x (m)     y (m)' )
2368601 FORMAT (2X,A1,1X,I2.2,6X,I2.2,5X,I5,5X,F8.2,2X,F8.2,5X,A)
[1]2369
2370 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.