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

Last change on this file since 3065 was 3065, checked in by Giersch, 6 years ago

New vertical stretching procedure has been introduced

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