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

Last change on this file since 1675 was 1675, checked in by gronemeier, 9 years ago

Bugfix: Definition of topography grid levels

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