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

Last change on this file since 1783 was 1783, checked in by raasch, 8 years ago

NetCDF routines modularized; new parameter netcdf_deflate; further changes in the pmc

  • Property svn:keywords set to Id
File size: 94.7 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!
[1691]16! Copyright 1997-2015 Leibniz Universitaet Hannover
[1036]17!--------------------------------------------------------------------------------!
18!
[254]19! Current revisions:
[1]20! -----------------
[1783]21! netcdf module and variable names changed, output of netcdf_deflate
[1698]22!
[1485]23! Former revisions:
24! -----------------
25! $Id: header.f90 1783 2016-03-06 18:36:17Z raasch $
26!
[1765]27! 1764 2016-02-28 12:45:19Z raasch
28! output of nesting informations
29!
[1698]30! 1697 2015-10-28 17:14:10Z raasch
31! small E- and F-FORMAT changes to avoid informative compiler messages about
32! insufficient field width
33!
[1692]34! 1691 2015-10-26 16:17:44Z maronga
35! Renamed prandtl_layer to constant_flux_layer, renames rif_min/rif_max to
36! zeta_min/zeta_max.
37!
[1683]38! 1682 2015-10-07 23:56:08Z knoop
39! Code annotations made doxygen readable
40!
[1676]41! 1675 2015-10-02 08:28:59Z gronemeier
42! Bugfix: Definition of topography grid levels
43!
[1662]44! 1660 2015-09-21 08:15:16Z gronemeier
45! Bugfix: Definition of building/street canyon height if vertical grid stretching
46!         starts below the maximum topography height.
47!
[1591]48! 1590 2015-05-08 13:56:27Z maronga
49! Bugfix: Added TRIM statements for character strings for LSM and radiation code
50!
[1586]51! 1585 2015-04-30 07:05:52Z maronga
52! Further output for radiation model(s).
53!
[1576]54! 1575 2015-03-27 09:56:27Z raasch
55! adjustments for psolver-queries, output of seed_follows_topography
56!
[1561]57! 1560 2015-03-06 10:48:54Z keck
58! output for recycling y shift
59!
[1558]60! 1557 2015-03-05 16:43:04Z suehring
61! output for monotonic limiter
62!
[1552]63! 1551 2015-03-03 14:18:16Z maronga
64! Added informal output for land surface model and radiation model. Removed typo.
65!
[1497]66! 1496 2014-12-02 17:25:50Z maronga
67! Renamed: "radiation -> "cloud_top_radiation"
68!
[1485]69! 1484 2014-10-21 10:53:05Z kanani
[1484]70! Changes due to new module structure of the plant canopy model:
71!   module plant_canopy_model_mod and output for new canopy model parameters
72!   (alpha_lad, beta_lad, lai_beta,...) added,
73!   drag_coefficient, leaf_surface_concentration and scalar_exchange_coefficient
74!   renamed to canopy_drag_coeff, leaf_surface_conc and leaf_scalar_exch_coeff,
75!   learde renamed leaf_area_density.
76! Bugfix: DO-WHILE-loop for lad header information additionally restricted
77! by maximum number of gradient levels (currently 10)
[1483]78!
79! 1482 2014-10-18 12:34:45Z raasch
80! information about calculated or predefined virtual processor topology adjusted
81!
[1469]82! 1468 2014-09-24 14:06:57Z maronga
83! Adapted for use on up to 6-digit processor cores
84!
[1430]85! 1429 2014-07-15 12:53:45Z knoop
86! header exended to provide ensemble_member_nr if specified
87!
[1377]88! 1376 2014-04-26 11:21:22Z boeske
89! Correction of typos
90!
[1366]91! 1365 2014-04-22 15:03:56Z boeske
92! New section 'Large scale forcing and nudging':
93! output of large scale forcing and nudging information,
94! new section for initial profiles created
95!
[1360]96! 1359 2014-04-11 17:15:14Z hoffmann
97! dt_sort_particles removed
98!
[1354]99! 1353 2014-04-08 15:21:23Z heinze
100! REAL constants provided with KIND-attribute
101!
[1329]102! 1327 2014-03-21 11:00:16Z raasch
103! parts concerning iso2d and avs output removed,
104! -netcdf output queries
105!
[1325]106! 1324 2014-03-21 09:13:16Z suehring
107! Bugfix: module spectrum added
108!
[1323]109! 1322 2014-03-20 16:38:49Z raasch
110! REAL functions provided with KIND-attribute,
111! some REAL constants defined as wp-kind
112!
[1321]113! 1320 2014-03-20 08:40:49Z raasch
[1320]114! ONLY-attribute added to USE-statements,
115! kind-parameters added to all INTEGER and REAL declaration statements,
116! kinds are defined in new module kinds,
117! revision history before 2012 removed,
118! comment fields (!:) to be used for variable explanations added to
119! all variable declaration statements
[1321]120!
[1309]121! 1308 2014-03-13 14:58:42Z fricke
122! output of the fixed number of output time levels
123! output_format adjusted for masked data if netcdf_data_format > 5
124!
[1300]125! 1299 2014-03-06 13:15:21Z heinze
126! output for using large_scale subsidence in combination
127! with large_scale_forcing
128! reformatting, more detailed explanations
129!
[1242]130! 1241 2013-10-30 11:36:58Z heinze
131! output for nudging + large scale forcing from external file
132!
[1217]133! 1216 2013-08-26 09:31:42Z raasch
134! output for transpose_compute_overlap
135!
[1213]136! 1212 2013-08-15 08:46:27Z raasch
137! output for poisfft_hybrid removed
138!
[1182]139! 1179 2013-06-14 05:57:58Z raasch
140! output of reference_state, use_reference renamed use_single_reference_value
141!
[1160]142! 1159 2013-05-21 11:58:22Z fricke
143! +use_cmax
144!
[1116]145! 1115 2013-03-26 18:16:16Z hoffmann
146! descriptions for Seifert-Beheng-cloud-physics-scheme added
147!
[1112]148! 1111 2013-03-08 23:54:10Z raasch
149! output of accelerator board information
150! ibc_p_b = 2 removed
151!
[1109]152! 1108 2013-03-05 07:03:32Z raasch
153! bugfix for r1106
154!
[1107]155! 1106 2013-03-04 05:31:38Z raasch
156! some format changes for coupled runs
157!
[1093]158! 1092 2013-02-02 11:24:22Z raasch
159! unused variables removed
160!
[1037]161! 1036 2012-10-22 13:43:42Z raasch
162! code put under GPL (PALM 3.9)
163!
[1035]164! 1031 2012-10-19 14:35:30Z raasch
165! output of netCDF data format modified
166!
[1017]167! 1015 2012-09-27 09:23:24Z raasch
[1365]168! output of Adjustment of mixing length to the Prandtl mixing length at first
[1017]169! grid point above ground removed
170!
[1004]171! 1003 2012-09-14 14:35:53Z raasch
172! output of information about equal/unequal subdomain size removed
173!
[1002]174! 1001 2012-09-13 14:08:46Z raasch
175! all actions concerning leapfrog- and upstream-spline-scheme removed
176!
[979]177! 978 2012-08-09 08:28:32Z fricke
178! -km_damp_max, outflow_damping_width
179! +pt_damping_factor, pt_damping_width
180! +z0h
181!
[965]182! 964 2012-07-26 09:14:24Z raasch
183! output of profil-related quantities removed
184!
[941]185! 940 2012-07-09 14:31:00Z raasch
186! Output in case of simulations for pure neutral stratification (no pt-equation
187! solved)
188!
[928]189! 927 2012-06-06 19:15:04Z raasch
190! output of masking_method for mg-solver
191!
[869]192! 868 2012-03-28 12:21:07Z raasch
193! translation velocity in Galilean transformation changed to 0.6 * ug
194!
[834]195! 833 2012-02-22 08:55:55Z maronga
196! Adjusted format for leaf area density
197!
[829]198! 828 2012-02-21 12:00:36Z raasch
199! output of dissipation_classes + radius_classes
200!
[826]201! 825 2012-02-19 03:03:44Z raasch
202! Output of cloud physics parameters/quantities complemented and restructured
203!
[1]204! Revision 1.1  1997/08/11 06:17:20  raasch
205! Initial revision
206!
207!
208! Description:
209! ------------
[1764]210!> Writing a header with all important information about the current run.
[1682]211!> This subroutine is called three times, two times at the beginning
212!> (writing information on files RUN_CONTROL and HEADER) and one time at the
213!> end of the run, then writing additional information about CPU-usage on file
214!> header.
[411]215!-----------------------------------------------------------------------------!
[1682]216 SUBROUTINE header
217 
[1]218
[1320]219    USE arrays_3d,                                                             &
[1660]220        ONLY:  pt_init, qsws, q_init, sa_init, shf, ug, vg, w_subs, zu, zw
[1320]221       
[1]222    USE control_parameters
[1320]223       
224    USE cloud_parameters,                                                      &
225        ONLY:  cp, curvature_solution_effects, c_sedimentation,                &
226               limiter_sedimentation, l_v, nc_const, r_d, ventilation_effect
227       
228    USE cpulog,                                                                &
229        ONLY:  log_point_s
230       
231    USE dvrp_variables,                                                        &
232        ONLY:  use_seperate_pe_for_dvrp_output
233       
234    USE grid_variables,                                                        &
235        ONLY:  dx, dy
236       
237    USE indices,                                                               &
238        ONLY:  mg_loc_ind, nnx, nny, nnz, nx, ny, nxl_mg, nxr_mg, nyn_mg,      &
239               nys_mg, nzt, nzt_mg
240       
241    USE kinds
[1551]242   
243    USE land_surface_model_mod,                                                &
244        ONLY:  conserve_water_content, dewfall, land_surface, nzb_soil,        &
245               nzt_soil, root_fraction, soil_moisture, soil_temperature,       &
246               soil_type, soil_type_name, veg_type, veg_type_name, zs
247 
[1320]248    USE model_1d,                                                              &
249        ONLY:  damp_level_ind_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
250       
[1783]251    USE netcdf_interface,                                                      &
252        ONLY:  netcdf_data_format, netcdf_data_format_string, netcdf_deflate
253
[1320]254    USE particle_attributes,                                                   &
255        ONLY:  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, collision_kernel,     &
256               density_ratio, dissipation_classes, dt_min_part, dt_prel,       &
[1359]257               dt_write_particle_data, end_time_prel,                          &
[1320]258               maximum_number_of_tailpoints, maximum_tailpoint_age,            &
259               minimum_tailpoint_distance, number_of_particle_groups,          &
260               particle_advection, particle_advection_start,                   &
261               particles_per_point, pdx, pdy, pdz,  psb, psl, psn, psr, pss,   &
262               pst, radius, radius_classes, random_start_position,             &
[1575]263               seed_follows_topography,                                        &
[1320]264               total_number_of_particles, use_particle_tails,                  &
265               use_sgs_for_particles, total_number_of_tails,                   &
266               vertical_particle_advection, write_particle_statistics
267       
[1]268    USE pegrid
[1484]269
270    USE plant_canopy_model_mod,                                                &
271        ONLY:  alpha_lad, beta_lad, calc_beta_lad_profile, canopy_drag_coeff,  &
272               canopy_mode, cthf, lad, lad_surface, lad_vertical_gradient,     &
273               lad_vertical_gradient_level, lad_vertical_gradient_level_ind,   &
274               lai_beta, leaf_scalar_exch_coeff, leaf_surface_conc, pch_index, &
275               plant_canopy
[1551]276
[1764]277    USE pmc_interface,                                                         &
278        ONLY:  cpl_id, cpl_parent_id, cpl_name, lower_left_coord_x,            &
279               lower_left_coord_y, nested_run, nesting_mode
280
[1551]281    USE radiation_model_mod,                                                   &
[1585]282        ONLY:  albedo, albedo_type, albedo_type_name, constant_albedo,         &
283               day_init, dt_radiation, lambda, lw_radiation, net_radiation,    &
284               radiation, radiation_scheme, sw_radiation, time_utc_init
[1324]285   
286    USE spectrum,                                                              &
287        ONLY:  comp_spectra_level, data_output_sp, plot_spectra_level,         &
288               spectra_direction
[1]289
290    IMPLICIT NONE
291
[1682]292    CHARACTER (LEN=1)  ::  prec                !<
[1320]293   
[1682]294    CHARACTER (LEN=2)  ::  do2d_mode           !<
[1320]295   
[1682]296    CHARACTER (LEN=5)  ::  section_chr         !<
[1320]297   
[1682]298    CHARACTER (LEN=10) ::  coor_chr            !<
299    CHARACTER (LEN=10) ::  host_chr            !<
[1320]300   
[1682]301    CHARACTER (LEN=16) ::  begin_chr           !<
[1320]302   
[1682]303    CHARACTER (LEN=26) ::  ver_rev             !<
[1320]304   
[1682]305    CHARACTER (LEN=40) ::  output_format       !<
[1320]306   
[1682]307    CHARACTER (LEN=70) ::  char1               !<
308    CHARACTER (LEN=70) ::  char2               !<
309    CHARACTER (LEN=70) ::  dopr_chr            !<
310    CHARACTER (LEN=70) ::  do2d_xy             !<
311    CHARACTER (LEN=70) ::  do2d_xz             !<
312    CHARACTER (LEN=70) ::  do2d_yz             !<
313    CHARACTER (LEN=70) ::  do3d_chr            !<
314    CHARACTER (LEN=70) ::  domask_chr          !<
315    CHARACTER (LEN=70) ::  run_classification  !<
[1320]316   
[1682]317    CHARACTER (LEN=85) ::  roben               !<
318    CHARACTER (LEN=85) ::  runten              !<
[1320]319   
[1682]320    CHARACTER (LEN=86) ::  coordinates         !<
321    CHARACTER (LEN=86) ::  gradients           !<
322    CHARACTER (LEN=86) ::  leaf_area_density   !<
323    CHARACTER (LEN=86) ::  roots               !<
324    CHARACTER (LEN=86) ::  slices              !<
325    CHARACTER (LEN=86) ::  temperatures        !<
326    CHARACTER (LEN=86) ::  ugcomponent         !<
327    CHARACTER (LEN=86) ::  vgcomponent         !<
[1]328
[1682]329    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)  !<
[410]330
[1682]331    INTEGER(iwp) ::  av        !<
332    INTEGER(iwp) ::  bh        !<
333    INTEGER(iwp) ::  blx       !<
334    INTEGER(iwp) ::  bly       !<
335    INTEGER(iwp) ::  bxl       !<
336    INTEGER(iwp) ::  bxr       !<
337    INTEGER(iwp) ::  byn       !<
338    INTEGER(iwp) ::  bys       !<
339    INTEGER(iwp) ::  ch        !<
340    INTEGER(iwp) ::  count     !<
341    INTEGER(iwp) ::  cwx       !<
342    INTEGER(iwp) ::  cwy       !<
343    INTEGER(iwp) ::  cxl       !<
344    INTEGER(iwp) ::  cxr       !<
345    INTEGER(iwp) ::  cyn       !<
346    INTEGER(iwp) ::  cys       !<
347    INTEGER(iwp) ::  dim       !<
348    INTEGER(iwp) ::  i         !<
349    INTEGER(iwp) ::  io        !<
350    INTEGER(iwp) ::  j         !<
351    INTEGER(iwp) ::  k         !<
352    INTEGER(iwp) ::  l         !<
353    INTEGER(iwp) ::  ll        !<
354    INTEGER(iwp) ::  mpi_type  !<
[1320]355   
[1682]356    REAL(wp) ::  canopy_height                    !< canopy height (in m)
357    REAL(wp) ::  cpuseconds_per_simulated_second  !<
[1]358
359!
360!-- Open the output file. At the end of the simulation, output is directed
361!-- to unit 19.
362    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
363         .NOT. simulated_time_at_begin /= simulated_time )  THEN
364       io = 15   !  header output on file RUN_CONTROL
365    ELSE
366       io = 19   !  header output on file HEADER
367    ENDIF
368    CALL check_open( io )
369
370!
371!-- At the end of the run, output file (HEADER) will be rewritten with
[1551]372!-- new information
[1]373    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
374
375!
376!-- Determine kind of model run
377    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
[1764]378       run_classification = 'restart run'
[328]379    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[1764]380       run_classification = 'run with cyclic fill of 3D - prerun data'
[147]381    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
[1764]382       run_classification = 'run without 1D - prerun'
[197]383    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
[1764]384       run_classification = 'run with 1D - prerun'
[197]385    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
[1764]386       run_classification = 'run initialized by user'
[1]387    ELSE
[254]388       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
389       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
[1]390    ENDIF
[1764]391    IF ( nested_run )  run_classification = 'nested ' // run_classification
[97]392    IF ( ocean )  THEN
393       run_classification = 'ocean - ' // run_classification
394    ELSE
395       run_classification = 'atmosphere - ' // run_classification
396    ENDIF
[1]397
398!
399!-- Run-identification, date, time, host
400    host_chr = host(1:10)
[75]401    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
[102]402    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
[291]403    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
404#if defined( __mpi2 )
405       mpi_type = 2
406#else
407       mpi_type = 1
408#endif
409       WRITE ( io, 101 )  mpi_type, coupling_mode
410    ENDIF
[1108]411#if defined( __parallel )
[1353]412    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]413       IF ( coupling_start_time > simulated_time_at_begin )  THEN
414          WRITE ( io, 109 )
415       ELSE
416          WRITE ( io, 114 )
417       ENDIF
418    ENDIF
[1108]419#endif
[1429]420    IF ( ensemble_member_nr /= 0 )  THEN
421       WRITE ( io, 512 )  run_date, run_identifier, run_time, runnr,           &
422                       ADJUSTR( host_chr ), ensemble_member_nr
423    ELSE
424       WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr,           &
[102]425                       ADJUSTR( host_chr )
[1429]426    ENDIF
[1]427#if defined( __parallel )
[1482]428    IF ( npex == -1  .AND.  npey == -1 )  THEN
[1]429       char1 = 'calculated'
430    ELSE
431       char1 = 'predefined'
432    ENDIF
433    IF ( threads_per_task == 1 )  THEN
[102]434       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
[1]435    ELSE
[102]436       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
[1]437                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
438    ENDIF
[1111]439    IF ( num_acc_per_node /= 0 )  WRITE ( io, 117 )  num_acc_per_node   
[1]440    IF ( ( host(1:3) == 'ibm'  .OR.  host(1:3) == 'nec'  .OR.    &
441           host(1:2) == 'lc'   .OR.  host(1:3) == 'dec' )  .AND. &
442         npex == -1  .AND.  pdims(2) == 1 )                      &
443    THEN
[102]444       WRITE ( io, 106 )
[1]445    ELSEIF ( pdims(2) == 1 )  THEN
[102]446       WRITE ( io, 107 )  'x'
[1]447    ELSEIF ( pdims(1) == 1 )  THEN
[102]448       WRITE ( io, 107 )  'y'
[1]449    ENDIF
[102]450    IF ( use_seperate_pe_for_dvrp_output )  WRITE ( io, 105 )
[759]451    IF ( numprocs /= maximum_parallel_io_streams )  THEN
452       WRITE ( io, 108 )  maximum_parallel_io_streams
453    ENDIF
[1111]454#else
455    IF ( num_acc_per_node /= 0 )  WRITE ( io, 120 )  num_acc_per_node
[1]456#endif
[1764]457
458!
459!-- Nesting informations
460    IF ( nested_run )  THEN
461       WRITE ( io, 600 )  cpl_id, TRIM( cpl_name ), cpl_parent_id,             &
462                          nesting_mode, lower_left_coord_x, lower_left_coord_y
463    ENDIF
[1]464    WRITE ( io, 99 )
465
466!
467!-- Numerical schemes
468    WRITE ( io, 110 )
469    IF ( psolver(1:7) == 'poisfft' )  THEN
470       WRITE ( io, 111 )  TRIM( fft_method )
[1216]471       IF ( transpose_compute_overlap )  WRITE( io, 115 )
[1]472    ELSEIF ( psolver == 'sor' )  THEN
473       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
[1575]474    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
475       WRITE ( io, 135 )  TRIM(psolver), cycle_mg, maximum_grid_level, ngsrb
[1]476       IF ( mg_cycles == -1 )  THEN
477          WRITE ( io, 140 )  residual_limit
478       ELSE
479          WRITE ( io, 141 )  mg_cycles
480       ENDIF
481       IF ( mg_switch_to_pe0_level == 0 )  THEN
482          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
483                             nzt_mg(1)
[197]484       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
[1]485          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
486                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
487                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
488                             nzt_mg(mg_switch_to_pe0_level),    &
489                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
490                             nzt_mg(1)
491       ENDIF
[927]492       IF ( masking_method )  WRITE ( io, 144 )
[1]493    ENDIF
494    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
495    THEN
496       WRITE ( io, 142 )
497    ENDIF
498
499    IF ( momentum_advec == 'pw-scheme' )  THEN
500       WRITE ( io, 113 )
[1299]501    ELSEIF (momentum_advec == 'ws-scheme' )  THEN
[667]502       WRITE ( io, 503 )
[1]503    ENDIF
504    IF ( scalar_advec == 'pw-scheme' )  THEN
505       WRITE ( io, 116 )
[667]506    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
507       WRITE ( io, 504 )
[1557]508    ELSEIF ( scalar_advec == 'ws-scheme-mono' )  THEN
509       WRITE ( io, 513 )
[1]510    ELSE
511       WRITE ( io, 118 )
512    ENDIF
[63]513
514    WRITE ( io, 139 )  TRIM( loop_optimization )
515
[1]516    IF ( galilei_transformation )  THEN
517       IF ( use_ug_for_galilei_tr )  THEN
[868]518          char1 = '0.6 * geostrophic wind'
[1]519       ELSE
520          char1 = 'mean wind in model domain'
521       ENDIF
522       IF ( simulated_time_at_begin == simulated_time )  THEN
523          char2 = 'at the start of the run'
524       ELSE
525          char2 = 'at the end of the run'
526       ENDIF
[1353]527       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ),                        &
528                          advected_distance_x/1000.0_wp,                       &
529                          advected_distance_y/1000.0_wp
[1]530    ENDIF
[1001]531    WRITE ( io, 122 )  timestep_scheme
[87]532    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
[1353]533    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
[108]534       IF ( .NOT. ocean )  THEN
535          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
536               rayleigh_damping_factor
537       ELSE
538          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
539               rayleigh_damping_factor
540       ENDIF
[1]541    ENDIF
[940]542    IF ( neutral )  WRITE ( io, 131 )  pt_surface
[75]543    IF ( humidity )  THEN
[1]544       IF ( .NOT. cloud_physics )  THEN
545          WRITE ( io, 129 )
546       ELSE
547          WRITE ( io, 130 )
548       ENDIF
549    ENDIF
550    IF ( passive_scalar )  WRITE ( io, 134 )
[240]551    IF ( conserve_volume_flow )  THEN
[241]552       WRITE ( io, 150 )  conserve_volume_flow_mode
553       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
554          WRITE ( io, 151 )  u_bulk, v_bulk
555       ENDIF
[240]556    ELSEIF ( dp_external )  THEN
557       IF ( dp_smooth )  THEN
[241]558          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
[240]559       ELSE
[241]560          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
[240]561       ENDIF
562    ENDIF
[1]563    WRITE ( io, 99 )
564
565!
[1551]566!-- Runtime and timestep information
[1]567    WRITE ( io, 200 )
568    IF ( .NOT. dt_fixed )  THEN
569       WRITE ( io, 201 )  dt_max, cfl_factor
570    ELSE
571       WRITE ( io, 202 )  dt
572    ENDIF
573    WRITE ( io, 203 )  simulated_time_at_begin, end_time
574
[1322]575    IF ( time_restart /= 9999999.9_wp  .AND. &
[1]576         simulated_time_at_begin == simulated_time )  THEN
[1322]577       IF ( dt_restart == 9999999.9_wp )  THEN
[1]578          WRITE ( io, 204 )  ' Restart at:       ',time_restart
579       ELSE
580          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
581       ENDIF
582    ENDIF
583
584    IF ( simulated_time_at_begin /= simulated_time )  THEN
585       i = MAX ( log_point_s(10)%counts, 1 )
[1353]586       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0_wp )  THEN
587          cpuseconds_per_simulated_second = 0.0_wp
[1]588       ELSE
589          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
590                                            ( simulated_time -    &
591                                              simulated_time_at_begin )
592       ENDIF
[1322]593       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum,      &
594                          log_point_s(10)%sum / REAL( i, KIND=wp ), &
[1]595                          cpuseconds_per_simulated_second
[1322]596       IF ( time_restart /= 9999999.9_wp  .AND.  time_restart < end_time )  THEN
597          IF ( dt_restart == 9999999.9_wp )  THEN
[1106]598             WRITE ( io, 204 )  ' Next restart at:     ',time_restart
[1]599          ELSE
[1106]600             WRITE ( io, 205 )  ' Next restart at:     ',time_restart, dt_restart
[1]601          ENDIF
602       ENDIF
603    ENDIF
604
[1324]605
[1]606!
[291]607!-- Start time for coupled runs, if independent precursor runs for atmosphere
[1106]608!-- and ocean are used or have been used. In this case, coupling_start_time
609!-- defines the time when the coupling is switched on.
[1353]610    IF ( coupling_start_time /= 0.0_wp )  THEN
[1106]611       WRITE ( io, 207 )  coupling_start_time
[291]612    ENDIF
613
614!
[1]615!-- Computational grid
[94]616    IF ( .NOT. ocean )  THEN
617       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
618       IF ( dz_stretch_level_index < nzt+1 )  THEN
619          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
620                             dz_stretch_factor, dz_max
621       ENDIF
622    ELSE
623       WRITE ( io, 250 )  dx, dy, dz, (nx+1)*dx, (ny+1)*dy, zu(0)
624       IF ( dz_stretch_level_index > 0 )  THEN
625          WRITE ( io, 252 )  dz_stretch_level, dz_stretch_level_index, &
626                             dz_stretch_factor, dz_max
627       ENDIF
[1]628    ENDIF
629    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ), &
630                       MIN( nnz+2, nzt+2 )
631    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
632
633!
[1365]634!-- Large scale forcing and nudging
635    WRITE ( io, 160 )
636    IF ( large_scale_forcing )  THEN
637       WRITE ( io, 162 )
638       WRITE ( io, 163 )
639
640       IF ( large_scale_subsidence )  THEN
641          IF ( .NOT. use_subsidence_tendencies )  THEN
642             WRITE ( io, 164 )
643          ELSE
644             WRITE ( io, 165 )
645          ENDIF
646       ENDIF
647
648       IF ( bc_pt_b == 'dirichlet' )  THEN
649          WRITE ( io, 180 )
650       ELSEIF ( bc_pt_b == 'neumann' )  THEN
651          WRITE ( io, 181 )
652       ENDIF
653
654       IF ( bc_q_b == 'dirichlet' )  THEN
655          WRITE ( io, 182 )
656       ELSEIF ( bc_q_b == 'neumann' )  THEN
657          WRITE ( io, 183 )
658       ENDIF
659
660       WRITE ( io, 167 )
661       IF ( nudging )  THEN
662          WRITE ( io, 170 )
663       ENDIF
664    ELSE
665       WRITE ( io, 161 )
666       WRITE ( io, 171 )
667    ENDIF
668    IF ( large_scale_subsidence )  THEN
669       WRITE ( io, 168 )
670       WRITE ( io, 169 )
671    ENDIF
672
673!
674!-- Profile for the large scale vertial velocity
675!-- Building output strings, starting with surface value
676    IF ( large_scale_subsidence )  THEN
677       temperatures = '   0.0'
678       gradients = '------'
679       slices = '     0'
680       coordinates = '   0.0'
681       i = 1
682       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
683
684          WRITE (coor_chr,'(E10.2,7X)')  &
685                                w_subs(subs_vertical_gradient_level_i(i))
686          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
687
688          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
689          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
690
691          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
692          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
693
694          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
695          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
696
697          IF ( i == 10 )  THEN
698             EXIT
699          ELSE
700             i = i + 1
701          ENDIF
702
703       ENDDO
704
705 
706       IF ( .NOT. large_scale_forcing )  THEN
707          WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
708                             TRIM( gradients ), TRIM( slices )
709       ENDIF
710
711
712    ENDIF
713
714!-- Profile of the geostrophic wind (component ug)
715!-- Building output strings
716    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
717    gradients = '------'
718    slices = '     0'
719    coordinates = '   0.0'
720    i = 1
721    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
722     
723       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
724       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
725
726       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
727       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
728
729       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
730       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
731
732       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
733       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
734
735       IF ( i == 10 )  THEN
736          EXIT
737       ELSE
738          i = i + 1
739       ENDIF
740
741    ENDDO
742
743    IF ( .NOT. large_scale_forcing )  THEN
744       WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
745                          TRIM( gradients ), TRIM( slices )
746    ENDIF
747
748!-- Profile of the geostrophic wind (component vg)
749!-- Building output strings
750    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
751    gradients = '------'
752    slices = '     0'
753    coordinates = '   0.0'
754    i = 1
755    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
756
757       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
758       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
759
760       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
761       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
762
763       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
764       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
765
766       WRITE (coor_chr,'(F6.1,1X)')  vg_vertical_gradient_level(i)
767       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
768
769       IF ( i == 10 )  THEN
770          EXIT
771       ELSE
772          i = i + 1
773       ENDIF
774 
775    ENDDO
776
777    IF ( .NOT. large_scale_forcing )  THEN
778       WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
779                          TRIM( gradients ), TRIM( slices )
780    ENDIF
781
782!
[1]783!-- Topography
784    WRITE ( io, 270 )  topography
785    SELECT CASE ( TRIM( topography ) )
786
787       CASE ( 'flat' )
788          ! no actions necessary
789
790       CASE ( 'single_building' )
791          blx = INT( building_length_x / dx )
792          bly = INT( building_length_y / dy )
[1675]793          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
794          IF ( ABS( zw(bh  ) - building_height ) == &
795               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
[1]796
[1322]797          IF ( building_wall_left == 9999999.9_wp )  THEN
[1]798             building_wall_left = ( nx + 1 - blx ) / 2 * dx
799          ENDIF
[1353]800          bxl = INT ( building_wall_left / dx + 0.5_wp )
[1]801          bxr = bxl + blx
802
[1322]803          IF ( building_wall_south == 9999999.9_wp )  THEN
[1]804             building_wall_south = ( ny + 1 - bly ) / 2 * dy
805          ENDIF
[1353]806          bys = INT ( building_wall_south / dy + 0.5_wp )
[1]807          byn = bys + bly
808
809          WRITE ( io, 271 )  building_length_x, building_length_y, &
810                             building_height, bxl, bxr, bys, byn
811
[240]812       CASE ( 'single_street_canyon' )
[1675]813          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
814          IF ( ABS( zw(ch  ) - canyon_height ) == &
815               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
[1322]816          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[240]817!
818!--          Street canyon in y direction
819             cwx = NINT( canyon_width_x / dx )
[1322]820             IF ( canyon_wall_left == 9999999.9_wp )  THEN
[240]821                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
822             ENDIF
823             cxl = NINT( canyon_wall_left / dx )
824             cxr = cxl + cwx
825             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
826
[1322]827          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[240]828!
829!--          Street canyon in x direction
830             cwy = NINT( canyon_width_y / dy )
[1322]831             IF ( canyon_wall_south == 9999999.9_wp )  THEN
[240]832                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
833             ENDIF
834             cys = NINT( canyon_wall_south / dy )
835             cyn = cys + cwy
836             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
837          ENDIF
838
[1]839    END SELECT
840
[256]841    IF ( TRIM( topography ) /= 'flat' )  THEN
842       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
843          IF ( TRIM( topography ) == 'single_building' .OR.  &
844               TRIM( topography ) == 'single_street_canyon' )  THEN
845             WRITE ( io, 278 )
846          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
847             WRITE ( io, 279 )
848          ENDIF
849       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
850          WRITE ( io, 278 )
851       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
852          WRITE ( io, 279 )
853       ENDIF
854    ENDIF
855
[1299]856    IF ( plant_canopy )  THEN
[1484]857   
858       canopy_height = pch_index * dz
[138]859
[1484]860       WRITE ( io, 280 )  canopy_mode, canopy_height, pch_index,               &
861                          canopy_drag_coeff
[1299]862       IF ( passive_scalar )  THEN
[1484]863          WRITE ( io, 281 )  leaf_scalar_exch_coeff,                           &
864                             leaf_surface_conc
[153]865       ENDIF
[138]866
[1]867!
[153]868!--    Heat flux at the top of vegetation
[1484]869       WRITE ( io, 282 )  cthf
[153]870
871!
[1484]872!--    Leaf area density profile, calculated either from given vertical
873!--    gradients or from beta probability density function.
874       IF (  .NOT.  calc_beta_lad_profile )  THEN
[138]875
[1484]876!--       Building output strings, starting with surface value
[1697]877          WRITE ( leaf_area_density, '(F7.4)' )  lad_surface
[1484]878          gradients = '------'
879          slices = '     0'
880          coordinates = '   0.0'
881          i = 1
882          DO  WHILE ( i < 11  .AND.  lad_vertical_gradient_level_ind(i) /= -9999 )
[138]883
[1484]884             WRITE (coor_chr,'(F7.2)')  lad(lad_vertical_gradient_level_ind(i))
885             leaf_area_density = TRIM( leaf_area_density ) // ' ' // TRIM( coor_chr )
886 
887             WRITE (coor_chr,'(F7.2)')  lad_vertical_gradient(i)
888             gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
[138]889
[1484]890             WRITE (coor_chr,'(I7)')  lad_vertical_gradient_level_ind(i)
891             slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
[138]892
[1484]893             WRITE (coor_chr,'(F7.1)')  lad_vertical_gradient_level(i)
894             coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
[138]895
[1484]896             i = i + 1
897          ENDDO
[138]898
[1484]899          WRITE ( io, 283 )  TRIM( coordinates ), TRIM( leaf_area_density ),              &
900                             TRIM( gradients ), TRIM( slices )
[138]901
[1484]902       ELSE
903       
[1697]904          WRITE ( leaf_area_density, '(F7.4)' )  lad_surface
[1484]905          coordinates = '   0.0'
906         
907          DO  k = 1, pch_index
908
909             WRITE (coor_chr,'(F7.2)')  lad(k)
910             leaf_area_density = TRIM( leaf_area_density ) // ' ' // TRIM( coor_chr )
911 
912             WRITE (coor_chr,'(F7.1)')  zu(k)
913             coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
914
915          ENDDO       
916
917          WRITE ( io, 284 ) TRIM( coordinates ), TRIM( leaf_area_density ), alpha_lad,    &
918                            beta_lad, lai_beta
919
920       ENDIF 
921
[138]922    ENDIF
923
[1484]924
[1551]925    IF ( land_surface )  THEN
926
927       temperatures = ''
928       gradients    = '' ! use for humidity here
929       coordinates  = '' ! use for height
930       roots        = '' ! use for root fraction
931       slices       = '' ! use for index
932
933       i = 1
934       DO i = nzb_soil, nzt_soil
935          WRITE (coor_chr,'(F10.2,7X)') soil_temperature(i)
936          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
937
938          WRITE (coor_chr,'(F10.2,7X)') soil_moisture(i)
939          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
940
941          WRITE (coor_chr,'(F10.2,7X)')  - zs(i)
942          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
943
944          WRITE (coor_chr,'(F10.2,7X)')  root_fraction(i)
945          roots = TRIM( roots ) // ' '  // TRIM( coor_chr )
946
947          WRITE (coor_chr,'(I10,7X)')  i
948          slices = TRIM( slices ) // ' '  // TRIM( coor_chr )
949
950
951       ENDDO
952
[138]953!
[1551]954!--    Write land surface model header
955       WRITE( io, 419 )
956       IF ( conserve_water_content )  THEN
957          WRITE( io, 440 )
958       ELSE
959          WRITE( io, 441 )
960       ENDIF
961
962       IF ( dewfall )  THEN
963          WRITE( io, 442 )
964       ELSE
965          WRITE( io, 443 )
966       ENDIF
967
[1590]968       WRITE( io, 438 ) TRIM( veg_type_name(veg_type) ),                       &
969                        TRIM (soil_type_name(soil_type) )
[1551]970       WRITE( io, 439 ) TRIM( coordinates ), TRIM( temperatures ),             &
971                        TRIM( gradients ), TRIM( roots ), TRIM( slices )
972
973
974    ENDIF
975
976    IF ( radiation )  THEN
977!
[1585]978!--    Write radiation model header
[1551]979       WRITE( io, 444 )
980
981       IF ( radiation_scheme == "constant" )  THEN
982          WRITE( io, 445 ) net_radiation
983       ELSEIF ( radiation_scheme == "clear-sky" )  THEN
984          WRITE( io, 446 )
[1585]985       ELSEIF ( radiation_scheme == "rrtmg" )  THEN
986          WRITE( io, 447 )
987          IF ( .NOT. lw_radiation )  WRITE( io, 458 )
988          IF ( .NOT. sw_radiation )  WRITE( io, 459 )
989       ENDIF
990
991       IF ( albedo_type == 0 )  THEN
992          WRITE( io, 448 ) albedo
[1551]993       ELSE
[1590]994          WRITE( io, 456 ) TRIM( albedo_type_name(albedo_type) )
[1551]995       ENDIF
[1585]996       IF ( constant_albedo )  THEN
997          WRITE( io, 457 )
998       ENDIF
[1551]999       WRITE( io, 449 ) dt_radiation
1000    ENDIF
1001
1002
1003!
[1]1004!-- Boundary conditions
1005    IF ( ibc_p_b == 0 )  THEN
1006       runten = 'p(0)     = 0      |'
1007    ELSEIF ( ibc_p_b == 1 )  THEN
1008       runten = 'p(0)     = p(1)   |'
1009    ENDIF
1010    IF ( ibc_p_t == 0 )  THEN
1011       roben  = 'p(nzt+1) = 0      |'
1012    ELSE
1013       roben  = 'p(nzt+1) = p(nzt) |'
1014    ENDIF
1015
1016    IF ( ibc_uv_b == 0 )  THEN
1017       runten = TRIM( runten ) // ' uv(0)     = -uv(1)                |'
1018    ELSE
1019       runten = TRIM( runten ) // ' uv(0)     = uv(1)                 |'
1020    ENDIF
[132]1021    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
1022       roben  = TRIM( roben  ) // ' uv(nzt+1) = 0                     |'
1023    ELSEIF ( ibc_uv_t == 0 )  THEN
[1]1024       roben  = TRIM( roben  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
1025    ELSE
1026       roben  = TRIM( roben  ) // ' uv(nzt+1) = uv(nzt)               |'
1027    ENDIF
1028
1029    IF ( ibc_pt_b == 0 )  THEN
[1551]1030       IF ( land_surface )  THEN
1031          runten = TRIM( runten ) // ' pt(0)     = from soil model'
1032       ELSE
1033          runten = TRIM( runten ) // ' pt(0)     = pt_surface'
1034       ENDIF
[102]1035    ELSEIF ( ibc_pt_b == 1 )  THEN
[1551]1036       runten = TRIM( runten ) // ' pt(0)     = pt(1)'
[102]1037    ELSEIF ( ibc_pt_b == 2 )  THEN
[1551]1038       runten = TRIM( runten ) // ' pt(0)     = from coupled model'
[1]1039    ENDIF
1040    IF ( ibc_pt_t == 0 )  THEN
[19]1041       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt_top'
1042    ELSEIF( ibc_pt_t == 1 )  THEN
1043       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt)'
1044    ELSEIF( ibc_pt_t == 2 )  THEN
1045       roben  = TRIM( roben  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
[667]1046
[1]1047    ENDIF
1048
1049    WRITE ( io, 300 )  runten, roben
1050
1051    IF ( .NOT. constant_diffusion )  THEN
1052       IF ( ibc_e_b == 1 )  THEN
1053          runten = 'e(0)     = e(1)'
1054       ELSE
1055          runten = 'e(0)     = e(1) = (u*/0.1)**2'
1056       ENDIF
1057       roben = 'e(nzt+1) = e(nzt) = e(nzt-1)'
1058
[97]1059       WRITE ( io, 301 )  'e', runten, roben       
[1]1060
1061    ENDIF
1062
[97]1063    IF ( ocean )  THEN
1064       runten = 'sa(0)    = sa(1)'
1065       IF ( ibc_sa_t == 0 )  THEN
1066          roben =  'sa(nzt+1) = sa_surface'
[1]1067       ELSE
[97]1068          roben =  'sa(nzt+1) = sa(nzt)'
[1]1069       ENDIF
[97]1070       WRITE ( io, 301 ) 'sa', runten, roben
1071    ENDIF
[1]1072
[97]1073    IF ( humidity )  THEN
1074       IF ( ibc_q_b == 0 )  THEN
[1551]1075          IF ( land_surface )  THEN
1076             runten = 'q(0)     = from soil model'
1077          ELSE
1078             runten = 'q(0)     = q_surface'
1079          ENDIF
1080
[97]1081       ELSE
1082          runten = 'q(0)     = q(1)'
1083       ENDIF
1084       IF ( ibc_q_t == 0 )  THEN
1085          roben =  'q(nzt)   = q_top'
1086       ELSE
1087          roben =  'q(nzt)   = q(nzt-1) + dq/dz'
1088       ENDIF
1089       WRITE ( io, 301 ) 'q', runten, roben
1090    ENDIF
[1]1091
[97]1092    IF ( passive_scalar )  THEN
1093       IF ( ibc_q_b == 0 )  THEN
1094          runten = 's(0)     = s_surface'
1095       ELSE
1096          runten = 's(0)     = s(1)'
1097       ENDIF
1098       IF ( ibc_q_t == 0 )  THEN
1099          roben =  's(nzt)   = s_top'
1100       ELSE
1101          roben =  's(nzt)   = s(nzt-1) + ds/dz'
1102       ENDIF
1103       WRITE ( io, 301 ) 's', runten, roben
[1]1104    ENDIF
1105
1106    IF ( use_surface_fluxes )  THEN
1107       WRITE ( io, 303 )
1108       IF ( constant_heatflux )  THEN
[1299]1109          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[1241]1110             WRITE ( io, 306 )  shf(0,0)
1111          ELSE
1112             WRITE ( io, 306 )  surface_heatflux
1113          ENDIF
[1]1114          IF ( random_heatflux )  WRITE ( io, 307 )
1115       ENDIF
[75]1116       IF ( humidity  .AND.  constant_waterflux )  THEN
[1299]1117          IF ( large_scale_forcing .AND. lsf_surf )  THEN
[1241]1118             WRITE ( io, 311 ) qsws(0,0)
1119          ELSE
1120             WRITE ( io, 311 ) surface_waterflux
1121          ENDIF
[1]1122       ENDIF
1123       IF ( passive_scalar  .AND.  constant_waterflux )  THEN
1124          WRITE ( io, 313 ) surface_waterflux
1125       ENDIF
1126    ENDIF
1127
[19]1128    IF ( use_top_fluxes )  THEN
1129       WRITE ( io, 304 )
[102]1130       IF ( coupling_mode == 'uncoupled' )  THEN
[151]1131          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
[102]1132          IF ( constant_top_heatflux )  THEN
1133             WRITE ( io, 306 )  top_heatflux
1134          ENDIF
1135       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1136          WRITE ( io, 316 )
[19]1137       ENDIF
[97]1138       IF ( ocean  .AND.  constant_top_salinityflux )  THEN
1139          WRITE ( io, 309 )  top_salinityflux
1140       ENDIF
[75]1141       IF ( humidity  .OR.  passive_scalar )  THEN
[19]1142          WRITE ( io, 315 )
1143       ENDIF
1144    ENDIF
1145
[1691]1146    IF ( constant_flux_layer )  THEN
1147       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length,                     &
1148                          z0h_factor*roughness_length, kappa,                  &
1149                          zeta_min, zeta_max
[1]1150       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
[75]1151       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
[1]1152          WRITE ( io, 312 )
1153       ENDIF
1154       IF ( passive_scalar  .AND.  .NOT. constant_waterflux )  THEN
1155          WRITE ( io, 314 )
1156       ENDIF
1157    ELSE
1158       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
[1691]1159          WRITE ( io, 310 )  zeta_min, zeta_max
[1]1160       ENDIF
1161    ENDIF
1162
1163    WRITE ( io, 317 )  bc_lr, bc_ns
[707]1164    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1159]1165       WRITE ( io, 318 )  use_cmax, pt_damping_width, pt_damping_factor       
[151]1166       IF ( turbulent_inflow )  THEN
[1560]1167          IF ( .NOT. recycling_yshift ) THEN
1168             WRITE ( io, 319 )  recycling_width, recycling_plane, &
1169                                inflow_damping_height, inflow_damping_width
1170          ELSE
1171             WRITE ( io, 322 )  recycling_width, recycling_plane, &
1172                                inflow_damping_height, inflow_damping_width
1173          END IF
[151]1174       ENDIF
[1]1175    ENDIF
1176
1177!
[1365]1178!-- Initial Profiles
1179    WRITE ( io, 321 )
1180!
1181!-- Initial wind profiles
1182    IF ( u_profile(1) /= 9999999.9_wp )  WRITE ( io, 427 )
1183
1184!
1185!-- Initial temperature profile
1186!-- Building output strings, starting with surface temperature
1187    WRITE ( temperatures, '(F6.2)' )  pt_surface
1188    gradients = '------'
1189    slices = '     0'
1190    coordinates = '   0.0'
1191    i = 1
1192    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
1193
1194       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
1195       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1196
1197       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
1198       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1199
1200       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
1201       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1202
1203       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
1204       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1205
1206       IF ( i == 10 )  THEN
1207          EXIT
1208       ELSE
1209          i = i + 1
1210       ENDIF
1211
1212    ENDDO
1213
1214    IF ( .NOT. nudging )  THEN
1215       WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1216                          TRIM( gradients ), TRIM( slices )
1217    ELSE
1218       WRITE ( io, 428 ) 
1219    ENDIF
1220
1221!
1222!-- Initial humidity profile
1223!-- Building output strings, starting with surface humidity
1224    IF ( humidity  .OR.  passive_scalar )  THEN
1225       WRITE ( temperatures, '(E8.1)' )  q_surface
1226       gradients = '--------'
1227       slices = '       0'
1228       coordinates = '     0.0'
1229       i = 1
1230       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1231         
1232          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1233          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1234
1235          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1236          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1237         
1238          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1239          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1240         
1241          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1242          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1243
1244          IF ( i == 10 )  THEN
1245             EXIT
1246          ELSE
1247             i = i + 1
1248          ENDIF
1249
1250       ENDDO
1251
1252       IF ( humidity )  THEN
1253          IF ( .NOT. nudging )  THEN
1254             WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ), &
1255                                TRIM( gradients ), TRIM( slices )
1256          ENDIF
1257       ELSE
1258          WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ), &
1259                             TRIM( gradients ), TRIM( slices )
1260       ENDIF
1261    ENDIF
1262
1263!
1264!-- Initial salinity profile
1265!-- Building output strings, starting with surface salinity
1266    IF ( ocean )  THEN
1267       WRITE ( temperatures, '(F6.2)' )  sa_surface
1268       gradients = '------'
1269       slices = '     0'
1270       coordinates = '   0.0'
1271       i = 1
1272       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1273
1274          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1275          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1276
1277          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1278          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1279
1280          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1281          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1282
1283          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1284          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1285
1286          IF ( i == 10 )  THEN
1287             EXIT
1288          ELSE
1289             i = i + 1
1290          ENDIF
1291
1292       ENDDO
1293
1294       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1295                          TRIM( gradients ), TRIM( slices )
1296    ENDIF
1297
1298
1299!
[1]1300!-- Listing of 1D-profiles
[151]1301    WRITE ( io, 325 )  dt_dopr_listing
[1353]1302    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1303       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1304    ENDIF
1305
1306!
1307!-- DATA output
1308    WRITE ( io, 330 )
[1353]1309    IF ( averaging_interval_pr /= 0.0_wp )  THEN
[151]1310       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
[1]1311    ENDIF
1312
1313!
1314!-- 1D-profiles
[346]1315    dopr_chr = 'Profile:'
[1]1316    IF ( dopr_n /= 0 )  THEN
1317       WRITE ( io, 331 )
1318
1319       output_format = ''
[1783]1320       output_format = netcdf_data_format_string
1321       IF ( netcdf_deflate == 0 )  THEN
1322          WRITE ( io, 344 )  output_format
1323       ELSE
1324          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1325       ENDIF
[1]1326
1327       DO  i = 1, dopr_n
1328          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
1329          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
1330             WRITE ( io, 332 )  dopr_chr
1331             dopr_chr = '       :'
1332          ENDIF
1333       ENDDO
1334
1335       IF ( dopr_chr /= '' )  THEN
1336          WRITE ( io, 332 )  dopr_chr
1337       ENDIF
1338       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
[1353]1339       IF ( skip_time_dopr /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dopr
[1]1340    ENDIF
1341
1342!
1343!-- 2D-arrays
1344    DO  av = 0, 1
1345
1346       i = 1
1347       do2d_xy = ''
1348       do2d_xz = ''
1349       do2d_yz = ''
1350       DO  WHILE ( do2d(av,i) /= ' ' )
1351
1352          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
1353          do2d_mode = do2d(av,i)(l-1:l)
1354
1355          SELECT CASE ( do2d_mode )
1356             CASE ( 'xy' )
1357                ll = LEN_TRIM( do2d_xy )
1358                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1359             CASE ( 'xz' )
1360                ll = LEN_TRIM( do2d_xz )
1361                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1362             CASE ( 'yz' )
1363                ll = LEN_TRIM( do2d_yz )
1364                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1365          END SELECT
1366
1367          i = i + 1
1368
1369       ENDDO
1370
1371       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
1372              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
[1327]1373              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) ) )  THEN
[1]1374
1375          IF (  av == 0 )  THEN
1376             WRITE ( io, 334 )  ''
1377          ELSE
1378             WRITE ( io, 334 )  '(time-averaged)'
1379          ENDIF
1380
1381          IF ( do2d_at_begin )  THEN
1382             begin_chr = 'and at the start'
1383          ELSE
1384             begin_chr = ''
1385          ENDIF
1386
1387          output_format = ''
[1783]1388          output_format = netcdf_data_format_string
1389          IF ( netcdf_deflate == 0 )  THEN
1390             WRITE ( io, 344 )  output_format
1391          ELSE
1392             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1393          ENDIF
[1]1394
1395          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
1396             i = 1
1397             slices = '/'
1398             coordinates = '/'
1399!
[1551]1400!--          Building strings with index and coordinate information of the
[1]1401!--          slices
1402             DO  WHILE ( section(i,1) /= -9999 )
1403
1404                WRITE (section_chr,'(I5)')  section(i,1)
1405                section_chr = ADJUSTL( section_chr )
1406                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1407
[206]1408                IF ( section(i,1) == -1 )  THEN
[1353]1409                   WRITE (coor_chr,'(F10.1)')  -1.0_wp
[206]1410                ELSE
1411                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
1412                ENDIF
[1]1413                coor_chr = ADJUSTL( coor_chr )
1414                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1415
1416                i = i + 1
1417             ENDDO
1418             IF ( av == 0 )  THEN
1419                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
1420                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
1421                                   TRIM( coordinates )
[1353]1422                IF ( skip_time_do2d_xy /= 0.0_wp )  THEN
[1]1423                   WRITE ( io, 339 )  skip_time_do2d_xy
1424                ENDIF
1425             ELSE
1426                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
1427                                   TRIM( begin_chr ), averaging_interval, &
1428                                   dt_averaging_input, 'k', TRIM( slices ), &
1429                                   TRIM( coordinates )
[1353]1430                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1431                   WRITE ( io, 339 )  skip_time_data_output_av
1432                ENDIF
1433             ENDIF
[1308]1434             IF ( netcdf_data_format > 4 )  THEN
1435                WRITE ( io, 352 )  ntdim_2d_xy(av)
1436             ELSE
1437                WRITE ( io, 353 )
1438             ENDIF
[1]1439          ENDIF
1440
1441          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
1442             i = 1
1443             slices = '/'
1444             coordinates = '/'
1445!
[1551]1446!--          Building strings with index and coordinate information of the
[1]1447!--          slices
1448             DO  WHILE ( section(i,2) /= -9999 )
1449
1450                WRITE (section_chr,'(I5)')  section(i,2)
1451                section_chr = ADJUSTL( section_chr )
1452                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1453
1454                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
1455                coor_chr = ADJUSTL( coor_chr )
1456                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1457
1458                i = i + 1
1459             ENDDO
1460             IF ( av == 0 )  THEN
1461                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
1462                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
1463                                   TRIM( coordinates )
[1353]1464                IF ( skip_time_do2d_xz /= 0.0_wp )  THEN
[1]1465                   WRITE ( io, 339 )  skip_time_do2d_xz
1466                ENDIF
1467             ELSE
1468                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
1469                                   TRIM( begin_chr ), averaging_interval, &
1470                                   dt_averaging_input, 'j', TRIM( slices ), &
1471                                   TRIM( coordinates )
[1353]1472                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1473                   WRITE ( io, 339 )  skip_time_data_output_av
1474                ENDIF
1475             ENDIF
[1308]1476             IF ( netcdf_data_format > 4 )  THEN
1477                WRITE ( io, 352 )  ntdim_2d_xz(av)
1478             ELSE
1479                WRITE ( io, 353 )
1480             ENDIF
[1]1481          ENDIF
1482
1483          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
1484             i = 1
1485             slices = '/'
1486             coordinates = '/'
1487!
[1551]1488!--          Building strings with index and coordinate information of the
[1]1489!--          slices
1490             DO  WHILE ( section(i,3) /= -9999 )
1491
1492                WRITE (section_chr,'(I5)')  section(i,3)
1493                section_chr = ADJUSTL( section_chr )
1494                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1495
1496                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
1497                coor_chr = ADJUSTL( coor_chr )
1498                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1499
1500                i = i + 1
1501             ENDDO
1502             IF ( av == 0 )  THEN
1503                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
1504                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
1505                                   TRIM( coordinates )
[1353]1506                IF ( skip_time_do2d_yz /= 0.0_wp )  THEN
[1]1507                   WRITE ( io, 339 )  skip_time_do2d_yz
1508                ENDIF
1509             ELSE
1510                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
1511                                   TRIM( begin_chr ), averaging_interval, &
1512                                   dt_averaging_input, 'i', TRIM( slices ), &
1513                                   TRIM( coordinates )
[1353]1514                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1515                   WRITE ( io, 339 )  skip_time_data_output_av
1516                ENDIF
1517             ENDIF
[1308]1518             IF ( netcdf_data_format > 4 )  THEN
1519                WRITE ( io, 352 )  ntdim_2d_yz(av)
1520             ELSE
1521                WRITE ( io, 353 )
1522             ENDIF
[1]1523          ENDIF
1524
1525       ENDIF
1526
1527    ENDDO
1528
1529!
1530!-- 3d-arrays
1531    DO  av = 0, 1
1532
1533       i = 1
1534       do3d_chr = ''
1535       DO  WHILE ( do3d(av,i) /= ' ' )
1536
1537          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
1538          i = i + 1
1539
1540       ENDDO
1541
1542       IF ( do3d_chr /= '' )  THEN
1543          IF ( av == 0 )  THEN
1544             WRITE ( io, 336 )  ''
1545          ELSE
1546             WRITE ( io, 336 )  '(time-averaged)'
1547          ENDIF
1548
[1783]1549          output_format = netcdf_data_format_string
1550          IF ( netcdf_deflate == 0 )  THEN
1551             WRITE ( io, 344 )  output_format
1552          ELSE
1553             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1554          ENDIF
[1]1555
1556          IF ( do3d_at_begin )  THEN
1557             begin_chr = 'and at the start'
1558          ELSE
1559             begin_chr = ''
1560          ENDIF
1561          IF ( av == 0 )  THEN
1562             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
1563                                zu(nz_do3d), nz_do3d
1564          ELSE
1565             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
1566                                TRIM( begin_chr ), averaging_interval, &
1567                                dt_averaging_input, zu(nz_do3d), nz_do3d
1568          ENDIF
1569
[1308]1570          IF ( netcdf_data_format > 4 )  THEN
1571             WRITE ( io, 352 )  ntdim_3d(av)
1572          ELSE
1573             WRITE ( io, 353 )
1574          ENDIF
1575
[1]1576          IF ( av == 0 )  THEN
[1353]1577             IF ( skip_time_do3d /= 0.0_wp )  THEN
[1]1578                WRITE ( io, 339 )  skip_time_do3d
1579             ENDIF
1580          ELSE
[1353]1581             IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[1]1582                WRITE ( io, 339 )  skip_time_data_output_av
1583             ENDIF
1584          ENDIF
1585
1586       ENDIF
1587
1588    ENDDO
1589
1590!
[410]1591!-- masked arrays
1592    IF ( masks > 0 )  WRITE ( io, 345 )  &
1593         mask_scale_x, mask_scale_y, mask_scale_z
1594    DO  mid = 1, masks
1595       DO  av = 0, 1
1596
1597          i = 1
1598          domask_chr = ''
1599          DO  WHILE ( domask(mid,av,i) /= ' ' )
1600             domask_chr = TRIM( domask_chr ) // ' ' //  &
1601                          TRIM( domask(mid,av,i) ) // ','
1602             i = i + 1
1603          ENDDO
1604
1605          IF ( domask_chr /= '' )  THEN
1606             IF ( av == 0 )  THEN
1607                WRITE ( io, 346 )  '', mid
1608             ELSE
1609                WRITE ( io, 346 )  ' (time-averaged)', mid
1610             ENDIF
1611
[1783]1612             output_format = netcdf_data_format_string
[1308]1613!--          Parallel output not implemented for mask data, hence
1614!--          output_format must be adjusted.
1615             IF ( netcdf_data_format == 5 ) output_format = 'netCDF4/HDF5'
1616             IF ( netcdf_data_format == 6 ) output_format = 'netCDF4/HDF5 classic'
[1783]1617             IF ( netcdf_deflate == 0 )  THEN
1618                WRITE ( io, 344 )  output_format
1619             ELSE
1620                WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1621             ENDIF
[410]1622
1623             IF ( av == 0 )  THEN
1624                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1625             ELSE
1626                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1627                                   averaging_interval, dt_averaging_input
1628             ENDIF
1629
1630             IF ( av == 0 )  THEN
[1353]1631                IF ( skip_time_domask(mid) /= 0.0_wp )  THEN
[410]1632                   WRITE ( io, 339 )  skip_time_domask(mid)
1633                ENDIF
1634             ELSE
[1353]1635                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
[410]1636                   WRITE ( io, 339 )  skip_time_data_output_av
1637                ENDIF
1638             ENDIF
1639!
1640!--          output locations
1641             DO  dim = 1, 3
[1353]1642                IF ( mask(mid,dim,1) >= 0.0_wp )  THEN
[410]1643                   count = 0
[1353]1644                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0_wp )
[410]1645                      count = count + 1
1646                   ENDDO
1647                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1648                                      mask(mid,dim,:count)
[1353]1649                ELSEIF ( mask_loop(mid,dim,1) < 0.0_wp .AND.  &
1650                         mask_loop(mid,dim,2) < 0.0_wp .AND.  &
1651                         mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1652                   WRITE ( io, 350 )  dir(dim), dir(dim)
[1353]1653                ELSEIF ( mask_loop(mid,dim,3) == 0.0_wp )  THEN
[410]1654                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1655                                      mask_loop(mid,dim,1:2)
1656                ELSE
1657                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1658                                      mask_loop(mid,dim,1:3)
1659                ENDIF
1660             ENDDO
1661          ENDIF
1662
1663       ENDDO
1664    ENDDO
1665
1666!
[1]1667!-- Timeseries
[1322]1668    IF ( dt_dots /= 9999999.9_wp )  THEN
[1]1669       WRITE ( io, 340 )
1670
[1783]1671       output_format = netcdf_data_format_string
1672       IF ( netcdf_deflate == 0 )  THEN
1673          WRITE ( io, 344 )  output_format
1674       ELSE
1675          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1676       ENDIF
[1]1677       WRITE ( io, 341 )  dt_dots
1678    ENDIF
1679
1680#if defined( __dvrp_graphics )
1681!
1682!-- Dvrp-output
[1322]1683    IF ( dt_dvrp /= 9999999.9_wp )  THEN
[1]1684       WRITE ( io, 360 )  dt_dvrp, TRIM( dvrp_output ), TRIM( dvrp_host ), &
1685                          TRIM( dvrp_username ), TRIM( dvrp_directory )
1686       i = 1
1687       l = 0
[336]1688       m = 0
[1]1689       DO WHILE ( mode_dvrp(i) /= ' ' )
1690          IF ( mode_dvrp(i)(1:10) == 'isosurface' )  THEN
[130]1691             READ ( mode_dvrp(i), '(10X,I2)' )  j
[1]1692             l = l + 1
1693             IF ( do3d(0,j) /= ' ' )  THEN
[336]1694                WRITE ( io, 361 )  TRIM( do3d(0,j) ), threshold(l), &
1695                                   isosurface_color(:,l)
[1]1696             ENDIF
1697          ELSEIF ( mode_dvrp(i)(1:6) == 'slicer' )  THEN
[130]1698             READ ( mode_dvrp(i), '(6X,I2)' )  j
[336]1699             m = m + 1
1700             IF ( do2d(0,j) /= ' ' )  THEN
1701                WRITE ( io, 362 )  TRIM( do2d(0,j) ), &
1702                                   slicer_range_limits_dvrp(:,m)
1703             ENDIF
[1]1704          ELSEIF ( mode_dvrp(i)(1:9) == 'particles' )  THEN
[336]1705             WRITE ( io, 363 )  dvrp_psize
1706             IF ( particle_dvrpsize /= 'none' )  THEN
1707                WRITE ( io, 364 )  'size', TRIM( particle_dvrpsize ), &
1708                                   dvrpsize_interval
1709             ENDIF
1710             IF ( particle_color /= 'none' )  THEN
1711                WRITE ( io, 364 )  'color', TRIM( particle_color ), &
1712                                   color_interval
1713             ENDIF
[1]1714          ENDIF
1715          i = i + 1
1716       ENDDO
[237]1717
[336]1718       WRITE ( io, 365 )  groundplate_color, superelevation_x, &
1719                          superelevation_y, superelevation, clip_dvrp_l, &
1720                          clip_dvrp_r, clip_dvrp_s, clip_dvrp_n
1721
1722       IF ( TRIM( topography ) /= 'flat' )  THEN
1723          WRITE ( io, 366 )  topography_color
1724          IF ( cluster_size > 1 )  THEN
1725             WRITE ( io, 367 )  cluster_size
1726          ENDIF
[237]1727       ENDIF
1728
[1]1729    ENDIF
1730#endif
1731
1732#if defined( __spectra )
1733!
1734!-- Spectra output
[1322]1735    IF ( dt_dosp /= 9999999.9_wp )  THEN
[1]1736       WRITE ( io, 370 )
1737
[1783]1738       output_format = netcdf_data_format_string
1739       IF ( netcdf_deflate == 0 )  THEN
1740          WRITE ( io, 344 )  output_format
1741       ELSE
1742          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1743       ENDIF
[1]1744       WRITE ( io, 371 )  dt_dosp
[1353]1745       IF ( skip_time_dosp /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dosp
[1]1746       WRITE ( io, 372 )  ( data_output_sp(i), i = 1,10 ),     &
1747                          ( spectra_direction(i), i = 1,10 ),  &
[189]1748                          ( comp_spectra_level(i), i = 1,100 ), &
1749                          ( plot_spectra_level(i), i = 1,100 ), &
[1]1750                          averaging_interval_sp, dt_averaging_input_pr
1751    ENDIF
1752#endif
1753
1754    WRITE ( io, 99 )
1755
1756!
1757!-- Physical quantities
1758    WRITE ( io, 400 )
1759
1760!
1761!-- Geostrophic parameters
[1551]1762    IF ( radiation .AND. radiation_scheme /= 'constant' )  THEN
1763       WRITE ( io, 417 )  lambda
1764    ENDIF
1765    WRITE ( io, 410 )  phi, omega, f, fs
[1]1766
1767!
1768!-- Other quantities
1769    WRITE ( io, 411 )  g
[1551]1770    IF ( radiation .AND. radiation_scheme /= 'constant' )  THEN
1771       WRITE ( io, 418 )  day_init, time_utc_init
1772    ENDIF
1773
[1179]1774    WRITE ( io, 412 )  TRIM( reference_state )
1775    IF ( use_single_reference_value )  THEN
[97]1776       IF ( ocean )  THEN
[1179]1777          WRITE ( io, 413 )  prho_reference
[97]1778       ELSE
[1179]1779          WRITE ( io, 414 )  pt_reference
[97]1780       ENDIF
1781    ENDIF
[1]1782
1783!
1784!-- Cloud physics parameters
[1299]1785    IF ( cloud_physics )  THEN
[57]1786       WRITE ( io, 415 )
1787       WRITE ( io, 416 ) surface_pressure, r_d, rho_surface, cp, l_v
[1115]1788       IF ( icloud_scheme == 0 )  THEN
[1353]1789          WRITE ( io, 510 ) 1.0E-6_wp * nc_const
[1115]1790          IF ( precipitation )  WRITE ( io, 511 ) c_sedimentation
1791       ENDIF
[1]1792    ENDIF
1793
1794!
[824]1795!-- Cloud physcis parameters / quantities / numerical methods
1796    WRITE ( io, 430 )
1797    IF ( humidity .AND. .NOT. cloud_physics .AND. .NOT. cloud_droplets)  THEN
1798       WRITE ( io, 431 )
1799    ELSEIF ( humidity  .AND.  cloud_physics )  THEN
1800       WRITE ( io, 432 )
[1496]1801       IF ( cloud_top_radiation )  WRITE ( io, 132 )
[1115]1802       IF ( icloud_scheme == 1 )  THEN
1803          IF ( precipitation )  WRITE ( io, 133 )
1804       ELSEIF ( icloud_scheme == 0 )  THEN
1805          IF ( drizzle )  WRITE ( io, 506 )
1806          IF ( precipitation )  THEN
1807             WRITE ( io, 505 )
1808             IF ( turbulence )  WRITE ( io, 507 )
1809             IF ( ventilation_effect )  WRITE ( io, 508 )
1810             IF ( limiter_sedimentation )  WRITE ( io, 509 )
1811          ENDIF
1812       ENDIF
[824]1813    ELSEIF ( humidity  .AND.  cloud_droplets )  THEN
1814       WRITE ( io, 433 )
1815       IF ( curvature_solution_effects )  WRITE ( io, 434 )
[825]1816       IF ( collision_kernel /= 'none' )  THEN
1817          WRITE ( io, 435 )  TRIM( collision_kernel )
[828]1818          IF ( collision_kernel(6:9) == 'fast' )  THEN
1819             WRITE ( io, 436 )  radius_classes, dissipation_classes
1820          ENDIF
[825]1821       ELSE
[828]1822          WRITE ( io, 437 )
[825]1823       ENDIF
[824]1824    ENDIF
1825
1826!
[1]1827!-- LES / turbulence parameters
1828    WRITE ( io, 450 )
1829
1830!--
1831! ... LES-constants used must still be added here
1832!--
1833    IF ( constant_diffusion )  THEN
1834       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1835                          prandtl_number
1836    ENDIF
1837    IF ( .NOT. constant_diffusion)  THEN
[1353]1838       IF ( e_init > 0.0_wp )  WRITE ( io, 455 )  e_init
1839       IF ( e_min > 0.0_wp )  WRITE ( io, 454 )  e_min
[1]1840       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1841    ENDIF
1842
1843!
1844!-- Special actions during the run
1845    WRITE ( io, 470 )
1846    IF ( create_disturbances )  THEN
1847       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1848                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1849                          zu(disturbance_level_ind_t), disturbance_level_ind_t
[707]1850       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1]1851          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1852       ELSE
1853          WRITE ( io, 473 )  disturbance_energy_limit
1854       ENDIF
1855       WRITE ( io, 474 )  TRIM( random_generator )
1856    ENDIF
[1353]1857    IF ( pt_surface_initial_change /= 0.0_wp )  THEN
[1]1858       WRITE ( io, 475 )  pt_surface_initial_change
1859    ENDIF
[1353]1860    IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1861       WRITE ( io, 476 )  q_surface_initial_change       
1862    ENDIF
[1353]1863    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
[1]1864       WRITE ( io, 477 )  q_surface_initial_change       
1865    ENDIF
1866
[60]1867    IF ( particle_advection )  THEN
[1]1868!
[60]1869!--    Particle attributes
1870       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
1871                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
[1359]1872                          end_time_prel
[60]1873       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
1874       IF ( random_start_position )  WRITE ( io, 481 )
[1575]1875       IF ( seed_follows_topography )  WRITE ( io, 496 )
[60]1876       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
1877       WRITE ( io, 495 )  total_number_of_particles
[824]1878       IF ( use_particle_tails  .AND.  maximum_number_of_tailpoints /= 0 )  THEN
[60]1879          WRITE ( io, 483 )  maximum_number_of_tailpoints
1880          IF ( minimum_tailpoint_distance /= 0 )  THEN
1881             WRITE ( io, 484 )  total_number_of_tails,      &
1882                                minimum_tailpoint_distance, &
1883                                maximum_tailpoint_age
1884          ENDIF
[1]1885       ENDIF
[1322]1886       IF ( dt_write_particle_data /= 9999999.9_wp )  THEN
[60]1887          WRITE ( io, 485 )  dt_write_particle_data
[1327]1888          IF ( netcdf_data_format > 1 )  THEN
1889             output_format = 'netcdf (64 bit offset) and binary'
[1]1890          ELSE
[1327]1891             output_format = 'netcdf and binary'
[1]1892          ENDIF
[1783]1893          IF ( netcdf_deflate == 0 )  THEN
1894             WRITE ( io, 344 )  output_format
1895          ELSE
1896             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1897          ENDIF
[1]1898       ENDIF
[1322]1899       IF ( dt_dopts /= 9999999.9_wp )  WRITE ( io, 494 )  dt_dopts
[60]1900       IF ( write_particle_statistics )  WRITE ( io, 486 )
[1]1901
[60]1902       WRITE ( io, 487 )  number_of_particle_groups
[1]1903
[60]1904       DO  i = 1, number_of_particle_groups
[1322]1905          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9_wp )  THEN
[1353]1906             WRITE ( io, 490 )  i, 0.0_wp
[60]1907             WRITE ( io, 492 )
[1]1908          ELSE
[60]1909             WRITE ( io, 490 )  i, radius(i)
[1353]1910             IF ( density_ratio(i) /= 0.0_wp )  THEN
[60]1911                WRITE ( io, 491 )  density_ratio(i)
1912             ELSE
1913                WRITE ( io, 492 )
1914             ENDIF
[1]1915          ENDIF
[60]1916          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
1917                             pdx(i), pdy(i), pdz(i)
[336]1918          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
[60]1919       ENDDO
[1]1920
[60]1921    ENDIF
[1]1922
[60]1923
[1]1924!
1925!-- Parameters of 1D-model
1926    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1927       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1928                          mixing_length_1d, dissipation_1d
1929       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1930          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1931       ENDIF
1932    ENDIF
1933
1934!
[1551]1935!-- User-defined information
[1]1936    CALL user_header( io )
1937
1938    WRITE ( io, 99 )
1939
1940!
1941!-- Write buffer contents to disc immediately
[82]1942    CALL local_flush( io )
[1]1943
1944!
1945!-- Here the FORMATs start
1946
1947 99 FORMAT (1X,78('-'))
[1468]1948100 FORMAT (/1X,'******************************',4X,44('-')/        &
1949            1X,'* ',A,' *',4X,A/                               &
1950            1X,'******************************',4X,44('-'))
1951101 FORMAT (35X,'coupled run using MPI-',I1,': ',A/ &
1952            35X,42('-'))
1953102 FORMAT (/' Date:                 ',A8,4X,'Run:       ',A20/      &
1954            ' Time:                 ',A8,4X,'Run-No.:   ',I2.2/     &
[1106]1955            ' Run on host:        ',A10)
[1]1956#if defined( __parallel )
[1468]1957103 FORMAT (' Number of PEs:',10X,I6,4X,'Processor grid (x,y): (',I4,',',I4, &
[1]1958              ')',1X,A)
[1468]1959104 FORMAT (' Number of PEs:',10X,I6,4X,'Tasks:',I4,'   threads per task:',I4/ &
1960              35X,'Processor grid (x,y): (',I4,',',I4,')',1X,A)
1961105 FORMAT (35X,'One additional PE is used to handle'/37X,'the dvrp output!')
1962106 FORMAT (35X,'A 1d-decomposition along x is forced'/ &
1963            35X,'because the job is running on an SMP-cluster')
1964107 FORMAT (35X,'A 1d-decomposition along ',A,' is used')
1965108 FORMAT (35X,'Max. # of parallel I/O streams is ',I5)
1966109 FORMAT (35X,'Precursor run for coupled atmos-ocean run'/ &
1967            35X,42('-'))
1968114 FORMAT (35X,'Coupled atmosphere-ocean run following'/ &
1969            35X,'independent precursor runs'/             &
1970            35X,42('-'))
[1111]1971117 FORMAT (' Accelerator boards / node:  ',I2)
[1]1972#endif
1973110 FORMAT (/' Numerical Schemes:'/ &
1974             ' -----------------'/)
1975111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1976112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
[1697]1977            '     Iterations (initial/other): ',I3,'/',I3,'  omega =',F6.3)
[1]1978113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1979                  ' or Upstream')
[1216]1980115 FORMAT ('     FFT and transpositions are overlapping')
[1]1981116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1982                  ' or Upstream')
1983118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
[1106]1984119 FORMAT (' --> Galilei-Transform applied to horizontal advection:'/ &
1985            '     translation velocity = ',A/ &
[1]1986            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
[1111]1987120 FORMAT (' Accelerator boards: ',8X,I2)
[1]1988122 FORMAT (' --> Time differencing scheme: ',A)
[108]1989123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
[1697]1990            '     maximum damping coefficient:',F6.3, ' 1/s')
[1]1991129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1992130 FORMAT (' --> Additional prognostic equation for the total water content')
[940]1993131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
1994                  F6.2, ' K assumed')
[824]1995132 FORMAT ('     Parameterization of long-wave radiation processes via'/ &
[1]1996            '     effective emissivity scheme')
[824]1997133 FORMAT ('     Precipitation parameterization via Kessler-Scheme')
[1]1998134 FORMAT (' --> Additional prognostic equation for a passive scalar')
[1575]1999135 FORMAT (' --> Solve perturbation pressure via ',A,' method (', &
[1]2000                  A,'-cycle)'/ &
2001            '     number of grid levels:                   ',I2/ &
2002            '     Gauss-Seidel red/black iterations:       ',I2)
2003136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
2004                  I3,')')
2005137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
2006            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
2007                  I3,')'/ &
2008            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
2009                  I3,')')
[63]2010139 FORMAT (' --> Loop optimization method: ',A)
[1]2011140 FORMAT ('     maximum residual allowed:                ',E10.3)
2012141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
2013142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
2014                  'step')
[87]2015143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
2016                  'kinetic energy')
[927]2017144 FORMAT ('     masking method is used')
[1]2018150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
[241]2019                  'conserved'/ &
2020            '     using the ',A,' mode')
2021151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
[306]2022152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
2023           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
2024           /'     starting from dp_level_b =', F8.3, 'm', A /)
[1365]2025160 FORMAT (//' Large scale forcing and nudging:'/ &
2026              ' -------------------------------'/)
2027161 FORMAT (' --> No large scale forcing from external is used (default) ')
2028162 FORMAT (' --> Large scale forcing from external file LSF_DATA is used: ')
2029163 FORMAT ('     - large scale advection tendencies ')
2030164 FORMAT ('     - large scale subsidence velocity w_subs ')
2031165 FORMAT ('     - large scale subsidence tendencies ')
2032167 FORMAT ('     - and geostrophic wind components ug and vg')
2033168 FORMAT (' --> Large-scale vertical motion is used in the ', &
[1299]2034                  'prognostic equation(s) for')
[1365]2035169 FORMAT ('     the scalar(s) only')
2036170 FORMAT (' --> Nudging is used')
2037171 FORMAT (' --> No nudging is used (default) ')
2038180 FORMAT ('     - prescribed surface values for temperature')
[1376]2039181 FORMAT ('     - prescribed surface fluxes for temperature')
2040182 FORMAT ('     - prescribed surface values for humidity')
[1365]2041183 FORMAT ('     - prescribed surface fluxes for humidity')
[1]2042200 FORMAT (//' Run time and time step information:'/ &
2043             ' ----------------------------------'/)
[1106]2044201 FORMAT ( ' Timestep:             variable     maximum value: ',F6.3,' s', &
[1697]2045             '    CFL-factor:',F5.2)
[1106]2046202 FORMAT ( ' Timestep:          dt = ',F6.3,' s'/)
2047203 FORMAT ( ' Start time:          ',F9.3,' s'/ &
2048             ' End time:            ',F9.3,' s')
[1]2049204 FORMAT ( A,F9.3,' s')
2050205 FORMAT ( A,F9.3,' s',5X,'restart every',17X,F9.3,' s')
[1106]2051206 FORMAT (/' Time reached:        ',F9.3,' s'/ &
2052             ' CPU-time used:       ',F9.3,' s     per timestep:               ', &
2053               '  ',F9.3,' s'/                                                    &
[1111]2054             '                                      per second of simulated tim', &
[1]2055               'e: ',F9.3,' s')
[1106]2056207 FORMAT ( ' Coupling start time: ',F9.3,' s')
[1]2057250 FORMAT (//' Computational grid and domain size:'/ &
2058              ' ----------------------------------'// &
2059              ' Grid length:      dx =    ',F7.3,' m    dy =    ',F7.3, &
2060              ' m    dz =    ',F7.3,' m'/ &
2061              ' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
2062              ' m  z(u) = ',F10.3,' m'/)
2063252 FORMAT (' dz constant up to ',F10.3,' m (k=',I4,'), then stretched by', &
[1697]2064              ' factor:',F6.3/ &
[1]2065            ' maximum dz not to be exceeded is dz_max = ',F10.3,' m'/)
2066254 FORMAT (' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
2067            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
2068260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
2069             ' degrees')
[1551]2070270 FORMAT (//' Topography information:'/ &
2071              ' ----------------------'// &
[1]2072              1X,'Topography: ',A)
2073271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
2074              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
2075                ' / ',I4)
[240]2076272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
2077              ' direction' / &
2078              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
2079              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
[256]2080278 FORMAT (' Topography grid definition convention:'/ &
2081            ' cell edge (staggered grid points'/  &
2082            ' (u in x-direction, v in y-direction))' /)
2083279 FORMAT (' Topography grid definition convention:'/ &
2084            ' cell center (scalar grid points)' /)
[138]2085280 FORMAT (//' Vegetation canopy (drag) model:'/ &
2086              ' ------------------------------'// &
2087              ' Canopy mode: ', A / &
[1484]2088              ' Canopy height: ',F6.2,'m (',I4,' grid points)' / &
[138]2089              ' Leaf drag coefficient: ',F6.2 /)
[1484]2090281 FORMAT (/ ' Scalar exchange coefficient: ',F6.2 / &
[153]2091              ' Scalar concentration at leaf surfaces in kg/m**3: ',F6.2 /)
2092282 FORMAT (' Predefined constant heatflux at the top of the vegetation: ',F6.2,' K m/s')
2093283 FORMAT (/ ' Characteristic levels of the leaf area density:'// &
[138]2094              ' Height:              ',A,'  m'/ &
2095              ' Leaf area density:   ',A,'  m**2/m**3'/ &
2096              ' Gradient:            ',A,'  m**2/m**4'/ &
2097              ' Gridpoint:           ',A)
[1484]2098284 FORMAT (//' Characteristic levels of the leaf area density and coefficients:'// &
2099              ' Height:              ',A,'  m'/ &
2100              ' Leaf area density:   ',A,'  m**2/m**3'/ &
2101              ' Coefficient alpha: ',F6.2 / &
2102              ' Coefficient beta: ',F6.2 / &
2103              ' Leaf area index: ',F6.2,'  m**2/m**2' /)
[138]2104               
[1]2105300 FORMAT (//' Boundary conditions:'/ &
2106             ' -------------------'// &
2107             '                     p                    uv             ', &
[1551]2108             '                     pt'// &
[1]2109             ' B. bound.: ',A/ &
2110             ' T. bound.: ',A)
[97]2111301 FORMAT (/'                     ',A// &
[1]2112             ' B. bound.: ',A/ &
2113             ' T. bound.: ',A)
[19]2114303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
2115304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
2116305 FORMAT (//'    Prandtl-Layer between bottom surface and first ', &
2117               'computational u,v-level:'// &
[1697]2118             '       zp = ',F6.2,' m   z0 =',F7.4,' m   z0h =',F8.5,&
2119             ' m   kappa =',F5.2/ &
2120             '       Rif value range:   ',F8.2,' <= rif <=',F6.2)
[97]2121306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
[1]2122307 FORMAT ('       Heatflux has a random normal distribution')
2123308 FORMAT ('       Predefined surface temperature')
[97]2124309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
[1]2125310 FORMAT (//'    1D-Model:'// &
2126             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
2127311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' m/s')
2128312 FORMAT ('       Predefined surface humidity')
2129313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
2130314 FORMAT ('       Predefined scalar value at the surface')
[19]2131315 FORMAT ('       Humidity / scalar flux at top surface is 0.0')
[102]2132316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
2133                    'atmosphere model')
[1]2134317 FORMAT (//' Lateral boundaries:'/ &
2135            '       left/right:  ',A/    &
2136            '       north/south: ',A)
[1159]2137318 FORMAT (/'       use_cmax: ',L1 / &
2138            '       pt damping layer width = ',F8.2,' m, pt ', &
[1697]2139                    'damping factor =',F7.4)
[151]2140319 FORMAT ('       turbulence recycling at inflow switched on'/ &
2141            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
2142            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
2143320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
[103]2144            '                                          v: ',F9.6,' m**2/s**2')
[1365]2145321 FORMAT (//' Initial profiles:'/ &
2146              ' ----------------')
[1560]2147322 FORMAT ('       turbulence recycling at inflow switched on'/ &
2148            '       y shift of the recycled inflow turbulence switched on'/ &
2149            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
[1592]2150            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m'/)
[151]2151325 FORMAT (//' List output:'/ &
[1]2152             ' -----------'//  &
2153            '    1D-Profiles:'/    &
2154            '       Output every             ',F8.2,' s')
[151]2155326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
[1]2156            '       Averaging input every    ',F8.2,' s')
2157330 FORMAT (//' Data output:'/ &
2158             ' -----------'/)
2159331 FORMAT (/'    1D-Profiles:')
2160332 FORMAT (/'       ',A)
2161333 FORMAT ('       Output every             ',F8.2,' s',/ &
2162            '       Time averaged over       ',F8.2,' s'/ &
2163            '       Averaging input every    ',F8.2,' s')
2164334 FORMAT (/'    2D-Arrays',A,':')
2165335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2166            '       Output every             ',F8.2,' s  ',A/ &
2167            '       Cross sections at ',A1,' = ',A/ &
2168            '       scalar-coordinates:   ',A,' m'/)
2169336 FORMAT (/'    3D-Arrays',A,':')
2170337 FORMAT (/'       Arrays: ',A/ &
2171            '       Output every             ',F8.2,' s  ',A/ &
2172            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
2173339 FORMAT ('       No output during initial ',F8.2,' s')
2174340 FORMAT (/'    Time series:')
2175341 FORMAT ('       Output every             ',F8.2,' s'/)
2176342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
2177            '       Output every             ',F8.2,' s  ',A/ &
2178            '       Time averaged over       ',F8.2,' s'/ &
2179            '       Averaging input every    ',F8.2,' s'/ &
2180            '       Cross sections at ',A1,' = ',A/ &
2181            '       scalar-coordinates:   ',A,' m'/)
2182343 FORMAT (/'       Arrays: ',A/ &
2183            '       Output every             ',F8.2,' s  ',A/ &
2184            '       Time averaged over       ',F8.2,' s'/ &
2185            '       Averaging input every    ',F8.2,' s'/ &
2186            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
[292]2187344 FORMAT ('       Output format: ',A/)
[410]2188345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
2189            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
2190            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
2191            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
2192346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
2193347 FORMAT ('       Variables: ',A/ &
2194            '       Output every             ',F8.2,' s')
2195348 FORMAT ('       Variables: ',A/ &
2196            '       Output every             ',F8.2,' s'/ &
2197            '       Time averaged over       ',F8.2,' s'/ &
2198            '       Averaging input every    ',F8.2,' s')
2199349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2200            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
2201            13('       ',8(F8.2,',')/) )
2202350 FORMAT (/'       Output locations in ',A,'-direction: ', &
2203            'all gridpoints along ',A,'-direction (default).' )
2204351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
2205            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
2206            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
[1313]2207352 FORMAT  (/'       Number of output time levels allowed: ',I3 /)
2208353 FORMAT  (/'       Number of output time levels allowed: unlimited' /)
[1783]2209354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
[1]2210#if defined( __dvrp_graphics )
2211360 FORMAT ('    Plot-Sequence with dvrp-software:'/ &
2212            '       Output every      ',F7.1,' s'/ &
2213            '       Output mode:      ',A/ &
2214            '       Host / User:      ',A,' / ',A/ &
2215            '       Directory:        ',A// &
2216            '       The sequence contains:')
[337]2217361 FORMAT (/'       Isosurface of "',A,'"    Threshold value: ', E12.3/ &
2218            '          Isosurface color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
2219362 FORMAT (/'       Slicer plane ',A/ &
[336]2220            '       Slicer limits: [',F6.2,',',F6.2,']')
[337]2221363 FORMAT (/'       Particles'/ &
[336]2222            '          particle size:  ',F7.2,' m')
2223364 FORMAT ('          particle ',A,' controlled by "',A,'" with interval [', &
2224                       F6.2,',',F6.2,']')
[337]2225365 FORMAT (/'       Groundplate color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)'/ &
[336]2226            '       Superelevation along (x,y,z): (',F4.1,',',F4.1,',',F4.1, &
2227                     ')'/ &
2228            '       Clipping limits: from x = ',F9.1,' m to x = ',F9.1,' m'/ &
2229            '                        from y = ',F9.1,' m to y = ',F9.1,' m')
[337]2230366 FORMAT (/'       Topography color: (',F4.2,',',F4.2,',',F4.2,') (R,G,B)')
[336]2231367 FORMAT ('       Polygon reduction for topography: cluster_size = ', I1)
[1]2232#endif
2233#if defined( __spectra )
2234370 FORMAT ('    Spectra:')
2235371 FORMAT ('       Output every ',F7.1,' s'/)
2236372 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
2237            '       Directions: ', 10(A5,',')/                         &
[189]2238            '       height levels  k = ', 20(I3,',')/                  &
2239            '                          ', 20(I3,',')/                  &
2240            '                          ', 20(I3,',')/                  &
2241            '                          ', 20(I3,',')/                  &
2242            '                          ', 19(I3,','),I3,'.'/           &
[1]2243            '       height levels selected for standard plot:'/        &
[189]2244            '                      k = ', 20(I3,',')/                  &
2245            '                          ', 20(I3,',')/                  &
2246            '                          ', 20(I3,',')/                  &
2247            '                          ', 20(I3,',')/                  &
2248            '                          ', 19(I3,','),I3,'.'/           &
[1]2249            '       Time averaged over ', F7.1, ' s,' /                &
2250            '       Profiles for the time averaging are taken every ', &
2251                    F6.1,' s')
2252#endif
2253400 FORMAT (//' Physical quantities:'/ &
2254              ' -------------------'/)
[1551]2255410 FORMAT ('    Geograph. latitude  :   phi    = ',F4.1,' degr'/   &
[1697]2256            '    Angular velocity    :   omega  =',E10.3,' rad/s'/  &
[1551]2257            '    Coriolis parameter  :   f      = ',F9.6,' 1/s'/    &
2258            '                            f*     = ',F9.6,' 1/s')
2259411 FORMAT (/'    Gravity             :   g      = ',F4.1,' m/s**2')
[1179]2260412 FORMAT (/'    Reference state used in buoyancy terms: ',A)
2261413 FORMAT ('       Reference density in buoyancy terms: ',F8.3,' kg/m**3')
2262414 FORMAT ('       Reference temperature in buoyancy terms: ',F8.4,' K')
[1551]2263415 FORMAT (/' Cloud physics parameters:'/ &
2264             ' ------------------------'/)
2265416 FORMAT ('    Surface pressure   :   p_0   = ',F7.2,' hPa'/      &
2266            '    Gas constant       :   R     = ',F5.1,' J/(kg K)'/ &
[1697]2267            '    Density of air     :   rho_0 =',F6.3,' kg/m**3'/  &
[1551]2268            '    Specific heat cap. :   c_p   = ',F6.1,' J/(kg K)'/ &
[1697]2269            '    Vapourization heat :   L_v   =',E9.2,' J/kg')
[1551]2270417 FORMAT ('    Geograph. longitude :   lambda = ',F4.1,' degr')
2271418 FORMAT (/'    Day of the year at model start :   day_init      =     ',I3 &
2272            /'    UTC time at model start        :   time_utc_init = ',F7.1' s')
2273419 FORMAT (//' Land surface model information:'/ &
2274              ' ------------------------------'/)
[1]2275420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
2276            '       Height:        ',A,'  m'/ &
2277            '       Temperature:   ',A,'  K'/ &
2278            '       Gradient:      ',A,'  K/100m'/ &
2279            '       Gridpoint:     ',A)
2280421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
2281            '       Height:      ',A,'  m'/ &
2282            '       Humidity:    ',A,'  kg/kg'/ &
2283            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
2284            '       Gridpoint:   ',A)
2285422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
2286            '       Height:                  ',A,'  m'/ &
2287            '       Scalar concentration:    ',A,'  kg/m**3'/ &
2288            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
2289            '       Gridpoint:               ',A)
2290423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
2291            '       Height:      ',A,'  m'/ &
2292            '       ug:          ',A,'  m/s'/ &
2293            '       Gradient:    ',A,'  1/100s'/ &
2294            '       Gridpoint:   ',A)
2295424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
2296            '       Height:      ',A,'  m'/ &
[97]2297            '       vg:          ',A,'  m/s'/ &
[1]2298            '       Gradient:    ',A,'  1/100s'/ &
2299            '       Gridpoint:   ',A)
[97]2300425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
2301            '       Height:     ',A,'  m'/ &
2302            '       Salinity:   ',A,'  psu'/ &
2303            '       Gradient:   ',A,'  psu/100m'/ &
2304            '       Gridpoint:  ',A)
[411]2305426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
2306            '       Height:      ',A,'  m'/ &
2307            '       w_subs:      ',A,'  m/s'/ &
2308            '       Gradient:    ',A,'  (m/s)/100m'/ &
2309            '       Gridpoint:   ',A)
[767]2310427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
2311                  ' profiles')
[1241]2312428 FORMAT (/'    Initial profiles (u, v, pt, q) are taken from file '/ &
2313             '    NUDGING_DATA')
[824]2314430 FORMAT (//' Cloud physics quantities / methods:'/ &
2315              ' ----------------------------------'/)
2316431 FORMAT ('    Humidity is treated as purely passive scalar (no condensati', &
2317                 'on)')
2318432 FORMAT ('    Bulk scheme with liquid water potential temperature and'/ &
2319            '    total water content is used.'/ &
2320            '    Condensation is parameterized via 0% - or 100% scheme.')
2321433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
2322                 'icle model')
2323434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
2324                 ' droplets < 1.0E-6 m')
[825]2325435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
[828]2326436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
2327                    'are used'/ &
2328            '          number of radius classes:       ',I3,'    interval ', &
2329                       '[1.0E-6,2.0E-4] m'/ &
2330            '          number of dissipation classes:   ',I2,'    interval ', &
2331                       '[0,1000] cm**2/s**3')
2332437 FORMAT ('    Droplet collision is switched off')
[1551]2333438 FORMAT (' --> Land surface type  : ',A,/ &
2334            ' --> Soil porosity type : ',A)
2335439 FORMAT (/'    Initial soil temperature and moisture profile:'// &
2336            '       Height:        ',A,'  m'/ &
2337            '       Temperature:   ',A,'  K'/ &
2338            '       Moisture:      ',A,'  m**3/m**3'/ &
2339            '       Root fraction: ',A,'  '/ &
2340            '       Gridpoint:     ',A)
2341440 FORMAT (/' --> Dewfall is allowed (default)')
2342441 FORMAT (' --> Dewfall is inhibited')
2343442 FORMAT (' --> Soil bottom is closed (water content is conserved, default)')
2344443 FORMAT (' --> Soil bottom is open (water content is not conserved)')
2345444 FORMAT (//' Radiation model information:'/                                 &
2346              ' ----------------------------'/)
2347445 FORMAT (' --> Using constant net radiation: net_radiation = ', F6.2, '  W/m**2')
2348446 FORMAT (' --> Simple radiation scheme for clear sky is used (no clouds,',  &
2349                   ' default)')
[1585]2350447 FORMAT (' --> RRTMG scheme is used')
[1697]2351448 FORMAT (/'     User-specific surface albedo: albedo =', F6.3)
[1585]2352449 FORMAT  ('     Timestep: dt_radiation = ', F5.2, '  s')
[1551]2353
[1]2354450 FORMAT (//' LES / Turbulence quantities:'/ &
2355              ' ---------------------------'/)
[824]2356451 FORMAT ('    Diffusion coefficients are constant:'/ &
2357            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
[1697]2358453 FORMAT ('    Mixing length is limited to',F5.2,' * z')
[824]2359454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
2360455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
[1585]2361456 FORMAT (/'    Albedo is set for land surface type: ', A)
2362457 FORMAT (/'    --> Albedo is fixed during the run')
2363458 FORMAT (/'    --> Longwave radiation is disabled')
2364459 FORMAT (/'    --> Shortwave radiation is disabled.')
[1]2365470 FORMAT (//' Actions during the simulation:'/ &
2366              ' -----------------------------'/)
[94]2367471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
[1697]2368            '    Disturbance amplitude           :    ',F5.2, ' m/s'/       &
[94]2369            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
2370            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
[1]2371472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
2372                 ' to i/j =',I4)
2373473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
[1697]2374                 F6.3, ' m**2/s**2')
[1]2375474 FORMAT ('    Random number generator used    : ',A/)
2376475 FORMAT ('    The surface temperature is increased (or decreased, ', &
2377                 'respectively, if'/ &
2378            '    the value is negative) by ',F5.2,' K at the beginning of the',&
2379                 ' 3D-simulation'/)
2380476 FORMAT ('    The surface humidity is increased (or decreased, ',&
2381                 'respectively, if the'/ &
2382            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
2383                 ' the 3D-simulation'/)
2384477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
2385                 'respectively, if the'/ &
2386            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
2387                 ' the 3D-simulation'/)
2388480 FORMAT ('    Particles:'/ &
2389            '    ---------'// &
2390            '       Particle advection is active (switched on at t = ', F7.1, &
2391                    ' s)'/ &
2392            '       Start of new particle generations every  ',F6.1,' s'/ &
2393            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
2394            '                            bottom:     ', A, ' top:         ', A/&
2395            '       Maximum particle age:                 ',F9.1,' s'/ &
[1359]2396            '       Advection stopped at t = ',F9.1,' s'/)
[1]2397481 FORMAT ('       Particles have random start positions'/)
[336]2398482 FORMAT ('          Particles are advected only horizontally'/)
[1]2399483 FORMAT ('       Particles have tails with a maximum of ',I3,' points')
2400484 FORMAT ('            Number of tails of the total domain: ',I10/ &
2401            '            Minimum distance between tailpoints: ',F8.2,' m'/ &
2402            '            Maximum age of the end of the tail:  ',F8.2,' s')
2403485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
2404486 FORMAT ('       Particle statistics are written on file'/)
2405487 FORMAT ('       Number of particle groups: ',I2/)
2406488 FORMAT ('       SGS velocity components are used for particle advection'/ &
[1697]2407            '          minimum timestep for advection:', F8.5/)
[1]2408489 FORMAT ('       Number of particles simultaneously released at each ', &
2409                    'point: ', I5/)
2410490 FORMAT ('       Particle group ',I2,':'/ &
2411            '          Particle radius: ',E10.3, 'm')
2412491 FORMAT ('          Particle inertia is activated'/ &
[1697]2413            '             density_ratio (rho_fluid/rho_particle) =',F6.3/)
[1]2414492 FORMAT ('          Particles are advected only passively (no inertia)'/)
2415493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
2416            '                                         y:',F8.1,' - ',F8.1,' m'/&
2417            '                                         z:',F8.1,' - ',F8.1,' m'/&
2418            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
2419                       ' m  dz = ',F8.1,' m'/)
2420494 FORMAT ('       Output of particle time series in NetCDF format every ', &
2421                    F8.2,' s'/)
2422495 FORMAT ('       Number of particles in total domain: ',I10/)
[1575]2423496 FORMAT ('       Initial vertical particle positions are interpreted ', &
2424                    'as relative to the given topography')
[1]2425500 FORMAT (//' 1D-Model parameters:'/                           &
2426              ' -------------------'//                           &
2427            '    Simulation time:                   ',F8.1,' s'/ &
2428            '    Run-controll output every:         ',F8.1,' s'/ &
2429            '    Vertical profile output every:     ',F8.1,' s'/ &
2430            '    Mixing length calculation:         ',A/         &
2431            '    Dissipation calculation:           ',A/)
2432502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
[667]2433503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
2434504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
[1115]2435505 FORMAT ('    Precipitation parameterization via Seifert-Beheng-Scheme')
2436506 FORMAT ('    Drizzle parameterization via Stokes law')
2437507 FORMAT ('    Turbulence effects on precipitation process')
2438508 FORMAT ('    Ventilation effects on evaporation of rain drops')
2439509 FORMAT ('    Slope limiter used for sedimentation process')
[1551]2440510 FORMAT ('    Droplet density    :   N_c   = ',F6.1,' 1/cm**3')
2441511 FORMAT ('    Sedimentation Courant number:                  '/&
[1697]2442            '                               C_s   =',F4.1,'        ')
[1429]2443512 FORMAT (/' Date:                 ',A8,6X,'Run:       ',A20/      &
2444            ' Time:                 ',A8,6X,'Run-No.:   ',I2.2/     &
2445            ' Run on host:        ',A10,6X,'En-No.:    ',I2.2)
[1557]2446513 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order ' // & 
2447            '+ monotonic adjustment')
[1764]2448600 FORMAT (/' Nesting informations:'/                                        &
2449            ' Nest id / name:                   ',I2.2,' / ',A,' (parent id: ',I2.2,')'/ &
2450            ' Nesting mode:                     ',A/ &
2451            ' Lower left corner coordinates:    ','x = ',F8.2,' m, y = ',F8.2,' m'/)
[1]2452
2453 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.