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

Last change on this file since 1630 was 1592, checked in by maronga, 10 years ago

removed typo in header

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