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

Last change on this file since 1832 was 1832, checked in by hoffmann, 8 years ago

last commit documented

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