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

Last change on this file since 4464 was 4444, checked in by raasch, 5 years ago

bugfix: cpp-directives for serial mode added

  • Property svn:keywords set to Id
File size: 76.9 KB
Line 
1! !> @file header.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
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!
17! Copyright 1997-2020 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: header.f90 4444 2020-03-05 15:59:50Z Giersch $
27! bugfix: cpp-directives for serial mode added
28!
29! 4360 2020-01-07 11:25:50Z suehring
30! Bugfix, character length too short, caused crash on NEC.
31!
32! 4309 2019-11-26 18:49:59Z suehring
33! replaced recycling_yshift by y_shift
34!
35! 4301 2019-11-22 12:09:09Z oliver.maas
36!
37! 4297 2019-11-21 10:37:50Z oliver.maas
38! Adjusted format for simulated time and related quantities
39!
40! 4297 2019-11-21 10:37:50Z oliver.maas
41! adjusted message to the changed parameter recycling_yshift
42!
43! 4227 2019-09-10 18:04:34Z gronemeier
44! implement new palm_date_time_mod
45!
46! 4223 2019-09-10 09:20:47Z gronemeier
47! Write information about rotation angle
48!
49! 4182 2019-08-22 15:20:23Z scharf
50! Corrected "Former revisions" section
51!
52! 4168 2019-08-16 13:50:17Z suehring
53! Replace function get_topography_top_index by topo_top_ind
54!
55! 4069 2019-07-01 14:05:51Z Giersch
56! Masked output running index mid has been introduced as a local variable to
57! avoid runtime error (Loop variable has been modified) in time_integration
58!
59! 4023 2019-06-12 13:20:01Z maronga
60! Renamed "coupling start time" to "spinup time"
61!
62! 4017 2019-06-06 12:16:46Z schwenkel
63! unused variable removed
64!
65! 3655 2019-01-07 16:51:22Z knoop
66! Implementation of the PALM module interface
67!
68! Revision 1.1  1997/08/11 06:17:20  raasch
69! Initial revision
70!
71!
72! Description:
73! ------------
74!> Writing a header with all important information about the current run.
75!> This subroutine is called three times, two times at the beginning
76!> (writing information on files RUN_CONTROL and HEADER) and one time at the
77!> end of the run, then writing additional information about CPU-usage on file
78!> header.
79!-----------------------------------------------------------------------------!
80 SUBROUTINE header
81 
82
83    USE arrays_3d,                                                             &
84        ONLY:  pt_init, q_init, s_init, sa_init, ug, vg, w_subs, zu, zw
85
86    USE basic_constants_and_equations_mod,                                     &
87        ONLY:  g, kappa
88
89    USE bulk_cloud_model_mod,                                                  &
90        ONLY:  bulk_cloud_model
91
92    USE control_parameters
93
94    USE cpulog,                                                                &
95        ONLY:  log_point_s
96
97    USE grid_variables,                                                        &
98        ONLY:  dx, dy
99
100    USE indices,                                                               &
101        ONLY:  mg_loc_ind, nnx, nny, nnz, nx, ny, nxl_mg, nxr_mg, nyn_mg,      &
102               nys_mg, nzt, nzt_mg, topo_top_ind
103
104    USE kinds
105
106    USE model_1d_mod,                                                          &
107        ONLY:  damp_level_ind_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
108
109    USE module_interface,                                                      &
110        ONLY:  module_interface_header
111
112    USE netcdf_interface,                                                      &
113        ONLY:  netcdf_data_format, netcdf_data_format_string, netcdf_deflate
114
115    USE ocean_mod,                                                             &
116        ONLY:  ibc_sa_t, prho_reference, sa_surface,                           &
117               sa_vertical_gradient, sa_vertical_gradient_level,               &
118               sa_vertical_gradient_level_ind
119
120    USE palm_date_time_mod,                                                    &
121        ONLY:  get_date_time
122
123    USE pegrid
124
125#if defined( __parallel )
126    USE pmc_handle_communicator,                                               &
127        ONLY:  pmc_get_model_info
128
129    USE pmc_interface,                                                         &
130        ONLY:  nested_run, nesting_datatransfer_mode, nesting_mode
131#endif
132
133    USE surface_mod,                                                           &
134        ONLY:  surf_def_h
135
136    USE turbulence_closure_mod,                                                &
137        ONLY:  rans_const_c, rans_const_sigma
138
139    IMPLICIT NONE
140
141   
142    CHARACTER (LEN=2)  ::  do2d_mode           !< mode of 2D data output (xy, xz, yz)
143   
144    CHARACTER (LEN=5)  ::  section_chr         !< string indicating grid information where to output 2D slices
145   
146    CHARACTER (LEN=10) ::  host_chr            !< string for hostname
147   
148    CHARACTER (LEN=16) ::  begin_chr           !< string indication start time for the data output
149    CHARACTER (LEN=16) ::  coor_chr            !< dummy string
150   
151    CHARACTER (LEN=26) ::  ver_rev             !< string for run identification
152
153#if defined( __parallel )
154    CHARACTER (LEN=32) ::  cpl_name            !< name of child domain (nesting mode only)
155#endif
156   
157    CHARACTER (LEN=40) ::  output_format       !< netcdf format
158       
159    CHARACTER (LEN=70) ::  char1               !< dummy varialbe used for various strings
160    CHARACTER (LEN=70) ::  char2               !< string containing informating about the advected distance in case of Galilei transformation
161    CHARACTER (LEN=23) ::  date_time_str       !< string for date and time information
162    CHARACTER (LEN=70) ::  dopr_chr            !< string indicating profile output variables
163    CHARACTER (LEN=70) ::  do2d_xy             !< string indicating 2D-xy output variables
164    CHARACTER (LEN=70) ::  do2d_xz             !< string indicating 2D-xz output variables
165    CHARACTER (LEN=70) ::  do2d_yz             !< string indicating 2D-yz output variables
166    CHARACTER (LEN=70) ::  do3d_chr            !< string indicating 3D output variables
167    CHARACTER (LEN=70) ::  domask_chr          !< string indicating masked output variables
168    CHARACTER (LEN=70) ::  run_classification  !< string classifying type of run, e.g. nested, coupled, etc.
169   
170    CHARACTER (LEN=85) ::  r_upper             !< string indicating model top boundary condition for various quantities
171    CHARACTER (LEN=85) ::  r_lower             !< string indicating bottom boundary condition for various quantities
172   
173    CHARACTER (LEN=86) ::  coordinates         !< string indicating height coordinates for profile-prescribed variables
174    CHARACTER (LEN=86) ::  gradients           !< string indicating gradients of profile-prescribed variables between the prescribed height coordinates
175    CHARACTER (LEN=86) ::  slices              !< string indicating grid coordinates of profile-prescribed subsidence velocity
176    CHARACTER (LEN=86) ::  temperatures        !< string indicating profile-prescribed subsidence velocities
177    CHARACTER (LEN=86) ::  ugcomponent         !< string indicating profile-prescribed geostrophic u-component
178    CHARACTER (LEN=86) ::  vgcomponent         !< string indicating profile-prescribed geostrophic v-component
179
180    CHARACTER (LEN=1), DIMENSION(1:3) ::  dir = (/ 'x', 'y', 'z' /)  !< string indicating masking steps along certain direction
181
182    INTEGER(iwp) ::  av             !< index indicating average output quantities
183    INTEGER(iwp) ::  bh             !< building height in generic single-building setup
184    INTEGER(iwp) ::  blx            !< building width in grid points along x in generic single-building setup
185    INTEGER(iwp) ::  bly            !< building width in grid points along y in generic single-building setup
186    INTEGER(iwp) ::  bxl            !< index for left building wall in generic single-building setup
187    INTEGER(iwp) ::  bxr            !< index for right building wall in generic single-building setup
188    INTEGER(iwp) ::  byn            !< index for north building wall in generic single-building setup
189    INTEGER(iwp) ::  bys            !< index for south building wall in generic single-building setup
190    INTEGER(iwp) ::  ch             !< canyon depth in generic street-canyon setup
191    INTEGER(iwp) ::  count          !< number of masked output locations
192    INTEGER(iwp) ::  cwx            !< canyon width along x in generic street-canyon setup
193    INTEGER(iwp) ::  cwy            !< canyon width along y in generic street-canyon setup
194    INTEGER(iwp) ::  cxl            !< index for left canyon wall in generic street-canyon setup
195    INTEGER(iwp) ::  cxr            !< index for right canyon wall in generic street-canyon setup
196    INTEGER(iwp) ::  cyn            !< index for north canyon wall in generic street-canyon setup
197    INTEGER(iwp) ::  cys            !< index for south canyon wall in generic street-canyon setup
198    INTEGER(iwp) ::  dim            !< running index for masking output locations
199    INTEGER(iwp) ::  i              !< running index for various loops
200    INTEGER(iwp) ::  io             !< file unit of HEADER file
201    INTEGER(iwp) ::  l              !< substring length
202    INTEGER(iwp) ::  ll             !< substring length
203    INTEGER(iwp) ::  mid            !< masked output running index
204#if defined( __parallel )
205    INTEGER(iwp) ::  cpl_parent_id  !< parent ID for the respective child model
206    INTEGER(iwp) ::  my_cpl_id      !< run id in a nested model setup
207    INTEGER(iwp) ::  n              !< running index over number of couplers in a nested model setup
208    INTEGER(iwp) ::  ncpl           !< number of coupler in a nested model setup
209    INTEGER(iwp) ::  npe_total      !< number of total PEs in a coupler (parent + child)
210#endif
211   
212
213    REAL(wp) ::  cpuseconds_per_simulated_second  !< CPU time (in s) per simulated second
214#if defined( __parallel )
215    REAL(wp) ::  lower_left_coord_x               !< x-coordinate of nest domain
216    REAL(wp) ::  lower_left_coord_y               !< y-coordinate of nest domain
217#endif
218
219!
220!-- Open the output file. At the end of the simulation, output is directed
221!-- to unit 19.
222    IF ( ( runnr == 0 .OR. force_print_header )  .AND. &
223         .NOT. simulated_time_at_begin /= simulated_time )  THEN
224       io = 15   !  header output on file RUN_CONTROL
225    ELSE
226       io = 19   !  header output on file HEADER
227    ENDIF
228    CALL check_open( io )
229
230!
231!-- At the end of the run, output file (HEADER) will be rewritten with
232!-- new information
233    IF ( io == 19 .AND. simulated_time_at_begin /= simulated_time ) REWIND( 19 )
234
235!
236!-- Determine kind of model run
237    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
238       run_classification = 'restart run'
239    ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
240       run_classification = 'run with cyclic fill of 3D - prerun data'
241    ELSEIF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 )  THEN
242       run_classification = 'run without 1D - prerun'
243    ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
244       run_classification = 'run with 1D - prerun'
245    ELSEIF ( INDEX( initializing_actions, 'inifor' ) /= 0 )  THEN
246       run_classification = 'run initialized with COSMO data'
247    ELSEIF ( INDEX( initializing_actions, 'by_user' ) /=0 )  THEN
248       run_classification = 'run initialized by user'
249    ELSEIF ( INDEX( initializing_actions, 'initialize_vortex' ) /=0 )  THEN
250       run_classification = 'run additionally initialized by a Rankine-vortex'
251    ELSEIF ( INDEX( initializing_actions, 'initialize_ptanom' ) /=0 )  THEN
252       run_classification = 'run additionally initialized by temperature anomaly'
253    ELSE
254       message_string = ' unknown action(s): ' // TRIM( initializing_actions )
255       CALL message( 'header', 'PA0191', 0, 0, 0, 6, 0 )
256    ENDIF
257#if defined( __parallel )
258    IF ( nested_run )  run_classification = 'nested ' // run_classification(1:63)
259#endif
260    IF ( ocean_mode )  THEN
261       run_classification = 'ocean - ' // run_classification(1:61)
262    ELSE
263       run_classification = 'atmosphere - ' // run_classification(1:57)
264    ENDIF
265
266!
267!-- Run-identification, date, time, host
268    host_chr = host(1:10)
269    ver_rev = TRIM( version ) // '  ' // TRIM( revision )
270    WRITE ( io, 100 )  ver_rev, TRIM( run_classification )
271    IF ( TRIM( coupling_mode ) /= 'uncoupled' )  THEN
272       WRITE ( io, 101 )  coupling_mode
273    ENDIF
274#if defined( __parallel )
275    IF ( coupling_start_time /= 0.0_wp  .AND. .NOT. spinup )  THEN
276       IF ( coupling_start_time > simulated_time_at_begin )  THEN
277          WRITE ( io, 109 )
278       ELSE
279          WRITE ( io, 114 )
280       ENDIF
281    ENDIF
282#endif
283    IF ( ensemble_member_nr /= 0 )  THEN
284       WRITE ( io, 512 )  run_date, run_identifier, run_time, runnr,           &
285                       ADJUSTR( host_chr ), ensemble_member_nr
286    ELSE
287       WRITE ( io, 102 )  run_date, run_identifier, run_time, runnr,           &
288                       ADJUSTR( host_chr )
289    ENDIF
290#if defined( __parallel )
291    IF ( npex == -1  .AND.  npey == -1 )  THEN
292       char1 = 'calculated'
293    ELSE
294       char1 = 'predefined'
295    ENDIF
296    IF ( threads_per_task == 1 )  THEN
297       WRITE ( io, 103 )  numprocs, pdims(1), pdims(2), TRIM( char1 )
298    ELSE
299       WRITE ( io, 104 )  numprocs*threads_per_task, numprocs, &
300                          threads_per_task, pdims(1), pdims(2), TRIM( char1 )
301    ENDIF
302
303    IF ( pdims(2) == 1 )  THEN
304       WRITE ( io, 107 )  'x'
305    ELSEIF ( pdims(1) == 1 )  THEN
306       WRITE ( io, 107 )  'y'
307    ENDIF
308    IF ( numprocs /= maximum_parallel_io_streams )  THEN
309       WRITE ( io, 108 )  maximum_parallel_io_streams
310    ENDIF
311#endif
312
313#if defined( __parallel )
314!
315!-- Nesting informations
316    IF ( nested_run )  THEN
317
318       WRITE ( io, 600 )  TRIM( nesting_mode ),                                &
319                          TRIM( nesting_datatransfer_mode )
320       CALL pmc_get_model_info( ncpl = ncpl, cpl_id = my_cpl_id )
321
322       DO  n = 1, ncpl
323          CALL pmc_get_model_info( request_for_cpl_id = n, cpl_name = cpl_name,&
324                                   cpl_parent_id = cpl_parent_id,              &
325                                   lower_left_x = lower_left_coord_x,          &
326                                   lower_left_y = lower_left_coord_y,          &
327                                   npe_total = npe_total )
328          IF ( n == my_cpl_id )  THEN
329             char1 = '*'
330          ELSE
331             char1 = ' '
332          ENDIF
333          WRITE ( io, 601 )  TRIM( char1 ), n, cpl_parent_id, npe_total,       &
334                             lower_left_coord_x, lower_left_coord_y,           &
335                             TRIM( cpl_name )
336       ENDDO
337
338    ENDIF
339#endif
340
341    WRITE ( io, 99 )
342
343!
344!-- Numerical schemes
345    WRITE ( io, 110 )
346    IF ( rans_mode )  THEN
347       WRITE ( io, 124 )  TRIM( turbulence_closure ), 'RANS'
348    ELSE
349       WRITE ( io, 124 )  TRIM( turbulence_closure ), 'LES'
350    ENDIF
351    WRITE ( io, 121 )  TRIM( approximation )
352    IF ( psolver(1:7) == 'poisfft' )  THEN
353       WRITE ( io, 111 )  TRIM( fft_method )
354       IF ( transpose_compute_overlap )  WRITE( io, 115 )
355    ELSEIF ( psolver == 'sor' )  THEN
356       WRITE ( io, 112 )  nsor_ini, nsor, omega_sor
357    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
358       WRITE ( io, 135 )  TRIM(psolver), cycle_mg, maximum_grid_level, ngsrb
359       IF ( mg_cycles == -1 )  THEN
360          WRITE ( io, 140 )  residual_limit
361       ELSE
362          WRITE ( io, 141 )  mg_cycles
363       ENDIF
364       IF ( mg_switch_to_pe0_level == 0 )  THEN
365          WRITE ( io, 136 )  nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
366                             nzt_mg(1)
367       ELSEIF (  mg_switch_to_pe0_level /= -1 )  THEN
368          WRITE ( io, 137 )  mg_switch_to_pe0_level,            &
369                             mg_loc_ind(2,0)-mg_loc_ind(1,0)+1, &
370                             mg_loc_ind(4,0)-mg_loc_ind(3,0)+1, &
371                             nzt_mg(mg_switch_to_pe0_level),    &
372                             nxr_mg(1)-nxl_mg(1)+1, nyn_mg(1)-nys_mg(1)+1, &
373                             nzt_mg(1)
374       ENDIF
375       IF ( psolver == 'multigrid_noopt' .AND. masking_method )  WRITE ( io, 144 )
376    ENDIF
377    IF ( call_psolver_at_all_substeps  .AND. timestep_scheme(1:5) == 'runge' ) &
378    THEN
379       WRITE ( io, 142 )
380    ENDIF
381
382    IF ( momentum_advec == 'pw-scheme' )  THEN
383       WRITE ( io, 113 )
384    ELSEIF (momentum_advec == 'ws-scheme' )  THEN
385       WRITE ( io, 503 )
386    ENDIF
387    IF ( scalar_advec == 'pw-scheme' )  THEN
388       WRITE ( io, 116 )
389    ELSEIF ( scalar_advec == 'ws-scheme' )  THEN
390       WRITE ( io, 504 )
391    ELSE
392       WRITE ( io, 118 )
393    ENDIF
394
395    WRITE ( io, 139 )  TRIM( loop_optimization )
396
397    IF ( galilei_transformation )  THEN
398       IF ( use_ug_for_galilei_tr )  THEN
399          char1 = '0.6 * geostrophic wind'
400       ELSE
401          char1 = 'mean wind in model domain'
402       ENDIF
403       IF ( simulated_time_at_begin == simulated_time )  THEN
404          char2 = 'at the start of the run'
405       ELSE
406          char2 = 'at the end of the run'
407       ENDIF
408       WRITE ( io, 119 )  TRIM( char1 ), TRIM( char2 ),                        &
409                          advected_distance_x/1000.0_wp,                       &
410                          advected_distance_y/1000.0_wp
411    ENDIF
412    WRITE ( io, 122 )  timestep_scheme
413    IF ( use_upstream_for_tke )  WRITE ( io, 143 )
414    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
415       IF ( .NOT. ocean_mode )  THEN
416          WRITE ( io, 123 )  'above', rayleigh_damping_height, &
417               rayleigh_damping_factor
418       ELSE
419          WRITE ( io, 123 )  'below', rayleigh_damping_height, &
420               rayleigh_damping_factor
421       ENDIF
422    ENDIF
423    IF ( neutral )  WRITE ( io, 131 )  pt_surface
424    IF ( humidity )  THEN
425       IF ( .NOT. bulk_cloud_model )  THEN
426          WRITE ( io, 129 )
427       ELSE
428          WRITE ( io, 130 )
429       ENDIF
430    ENDIF
431    IF ( passive_scalar )  WRITE ( io, 134 )
432    IF ( conserve_volume_flow )  THEN
433       WRITE ( io, 150 )  conserve_volume_flow_mode
434       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
435          WRITE ( io, 151 )  u_bulk, v_bulk
436       ENDIF
437    ELSEIF ( dp_external )  THEN
438       IF ( dp_smooth )  THEN
439          WRITE ( io, 152 )  dpdxy, dp_level_b, ', vertically smoothed.'
440       ELSE
441          WRITE ( io, 152 )  dpdxy, dp_level_b, '.'
442       ENDIF
443    ENDIF
444    WRITE ( io, 99 )
445
446!
447!-- Runtime and timestep information
448    WRITE ( io, 200 )
449    IF ( .NOT. dt_fixed )  THEN
450       WRITE ( io, 201 )  dt_max, cfl_factor
451    ELSE
452       WRITE ( io, 202 )  dt
453    ENDIF
454    WRITE ( io, 203 )  simulated_time_at_begin, end_time
455
456    IF ( time_restart /= 9999999.9_wp  .AND. &
457         simulated_time_at_begin == simulated_time )  THEN
458       IF ( dt_restart == 9999999.9_wp )  THEN
459          WRITE ( io, 204 )  ' Restart at:       ',time_restart
460       ELSE
461          WRITE ( io, 205 )  ' Restart at:       ',time_restart, dt_restart
462       ENDIF
463    ENDIF
464
465    IF ( simulated_time_at_begin /= simulated_time )  THEN
466       i = MAX ( log_point_s(10)%counts, 1 )
467       IF ( ( simulated_time - simulated_time_at_begin ) == 0.0_wp )  THEN
468          cpuseconds_per_simulated_second = 0.0_wp
469       ELSE
470          cpuseconds_per_simulated_second = log_point_s(10)%sum / &
471                                            ( simulated_time -    &
472                                              simulated_time_at_begin )
473       ENDIF
474       WRITE ( io, 206 )  simulated_time, log_point_s(10)%sum,      &
475                          log_point_s(10)%sum / REAL( i, KIND=wp ), &
476                          cpuseconds_per_simulated_second
477       IF ( time_restart /= 9999999.9_wp  .AND.  time_restart < end_time )  THEN
478          IF ( dt_restart == 9999999.9_wp )  THEN
479             WRITE ( io, 204 )  ' Next restart at:     ',time_restart
480          ELSE
481             WRITE ( io, 205 )  ' Next restart at:     ',time_restart, dt_restart
482          ENDIF
483       ENDIF
484    ENDIF
485
486
487!
488!-- Start time for coupled runs, if independent precursor runs for atmosphere
489!-- and ocean are used or have been used. In this case, coupling_start_time
490!-- defines the time when the coupling is switched on.
491    IF ( coupling_start_time /= 0.0_wp )  THEN
492       WRITE ( io, 207 )  coupling_start_time
493    ENDIF
494
495!
496!-- Computational grid
497    IF ( .NOT. ocean_mode )  THEN
498       WRITE ( io, 250 )  dx, dy
499       
500       DO i = 1, number_stretch_level_start+1
501          WRITE ( io, 253 )  i, dz(i)
502       ENDDO
503       
504       WRITE( io, 251 ) (nx+1)*dx, (ny+1)*dy, zu(nzt+1)
505       
506       IF ( ANY( dz_stretch_level_start_index < nzt+1 ) )  THEN
507          WRITE( io, '(A)', advance='no') ' Vertical stretching starts at height:'
508          DO i = 1, number_stretch_level_start
509             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_start(i), ' m,'
510          ENDDO
511          WRITE( io, '(/,A)', advance='no') ' Vertical stretching starts at index: '
512          DO i = 1, number_stretch_level_start
513             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_start_index(i), ','
514          ENDDO
515          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at height:  '
516          DO i = 1, number_stretch_level_start
517             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_end(i), ' m,'
518          ENDDO
519          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at index:   '
520          DO i = 1, number_stretch_level_start
521             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_end_index(i), ','
522          ENDDO
523          WRITE( io, '(/,A)', advance='no') ' Factor used for stretching:          '
524          DO i = 1, number_stretch_level_start
525             WRITE ( io, '(F12.3,A1)', advance='no' )  dz_stretch_factor_array(i), ','
526          ENDDO
527       ENDIF
528       
529    ELSE
530       WRITE ( io, 250 )  dx, dy
531       DO i = 1, number_stretch_level_start+1
532          WRITE ( io, 253 )  i, dz(i)
533       ENDDO
534       
535       WRITE ( io, 251 ) (nx+1)*dx, (ny+1)*dy, zu(0)
536       
537       IF ( ANY( dz_stretch_level_start_index > 0 ) )  THEN
538          WRITE( io, '(A)', advance='no') ' Vertical stretching starts at height:'
539          DO i = 1, number_stretch_level_start
540             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_start(i), ' m,'
541          ENDDO
542          WRITE( io, '(/,A)', advance='no') ' Vertical stretching starts at index: '
543          DO i = 1, number_stretch_level_start
544             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_start_index(i), ','
545          ENDDO
546          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at height:  '
547          DO i = 1, number_stretch_level_start
548             WRITE ( io, '(F10.1,A3)', advance='no' )  dz_stretch_level_end(i), ' m,'
549          ENDDO
550          WRITE( io, '(/,A)', advance='no') ' Vertical stretching ends at index:   '
551          DO i = 1, number_stretch_level_start
552             WRITE ( io, '(I12,A1)', advance='no' )  dz_stretch_level_end_index(i), ','
553          ENDDO
554          WRITE( io, '(/,A)', advance='no') ' Factor used for stretching:          '
555          DO i = 1, number_stretch_level_start
556             WRITE ( io, '(F12.3,A1)', advance='no' )  dz_stretch_factor_array(i), ','
557          ENDDO
558       ENDIF
559    ENDIF
560    WRITE ( io, 254 )  nx, ny, nzt+1, MIN( nnx, nx+1 ), MIN( nny, ny+1 ),      &
561                       MIN( nnz+2, nzt+2 )
562    IF ( sloping_surface )  WRITE ( io, 260 )  alpha_surface
563
564!
565!-- Profile for the large scale vertial velocity
566!-- Building output strings, starting with surface value
567    IF ( large_scale_subsidence )  THEN
568       temperatures = '   0.0'
569       gradients = '------'
570       slices = '     0'
571       coordinates = '   0.0'
572       i = 1
573       DO  WHILE ( subs_vertical_gradient_level_i(i) /= -9999 )
574
575          WRITE (coor_chr,'(E10.2,7X)')  &
576                                w_subs(subs_vertical_gradient_level_i(i))
577          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
578
579          WRITE (coor_chr,'(E10.2,7X)')  subs_vertical_gradient(i)
580          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
581
582          WRITE (coor_chr,'(I10,7X)')  subs_vertical_gradient_level_i(i)
583          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
584
585          WRITE (coor_chr,'(F10.2,7X)')  subs_vertical_gradient_level(i)
586          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
587
588          IF ( i == 10 )  THEN
589             EXIT
590          ELSE
591             i = i + 1
592          ENDIF
593
594       ENDDO
595
596 
597       IF ( .NOT. large_scale_forcing )  THEN
598          WRITE ( io, 426 )  TRIM( coordinates ), TRIM( temperatures ), &
599                             TRIM( gradients ), TRIM( slices )
600       ENDIF
601
602
603    ENDIF
604
605!-- Profile of the geostrophic wind (component ug)
606!-- Building output strings
607    WRITE ( ugcomponent, '(F6.2)' )  ug_surface
608    gradients = '------'
609    slices = '     0'
610    coordinates = '   0.0'
611    i = 1
612    DO  WHILE ( ug_vertical_gradient_level_ind(i) /= -9999 )
613     
614       WRITE (coor_chr,'(F6.2,1X)')  ug(ug_vertical_gradient_level_ind(i))
615       ugcomponent = TRIM( ugcomponent ) // '  ' // TRIM( coor_chr )
616
617       WRITE (coor_chr,'(F6.2,1X)')  ug_vertical_gradient(i)
618       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
619
620       WRITE (coor_chr,'(I6,1X)')  ug_vertical_gradient_level_ind(i)
621       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
622
623       WRITE (coor_chr,'(F6.1,1X)')  ug_vertical_gradient_level(i)
624       coordinates = TRIM( coordinates ) // '  ' // TRIM( coor_chr )
625
626       IF ( i == 10 )  THEN
627          EXIT
628       ELSE
629          i = i + 1
630       ENDIF
631
632    ENDDO
633
634    IF ( .NOT. large_scale_forcing )  THEN
635       WRITE ( io, 423 )  TRIM( coordinates ), TRIM( ugcomponent ), &
636                          TRIM( gradients ), TRIM( slices )
637    ENDIF
638
639!-- Profile of the geostrophic wind (component vg)
640!-- Building output strings
641    WRITE ( vgcomponent, '(F6.2)' )  vg_surface
642    gradients = '------'
643    slices = '     0'
644    coordinates = '   0.0'
645    i = 1
646    DO  WHILE ( vg_vertical_gradient_level_ind(i) /= -9999 )
647
648       WRITE (coor_chr,'(F6.2,1X)')  vg(vg_vertical_gradient_level_ind(i))
649       vgcomponent = TRIM( vgcomponent ) // '  ' // TRIM( coor_chr )
650
651       WRITE (coor_chr,'(F6.2,1X)')  vg_vertical_gradient(i)
652       gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
653
654       WRITE (coor_chr,'(I6,1X)')  vg_vertical_gradient_level_ind(i)
655       slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
656
657       WRITE (coor_chr,'(F6.1,1X)')  vg_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    IF ( .NOT. large_scale_forcing )  THEN
669       WRITE ( io, 424 )  TRIM( coordinates ), TRIM( vgcomponent ), &
670                          TRIM( gradients ), TRIM( slices )
671    ENDIF
672
673!
674!-- Topography
675    WRITE ( io, 270 )  topography
676    SELECT CASE ( TRIM( topography ) )
677
678       CASE ( 'flat' )
679          ! no actions necessary
680
681       CASE ( 'single_building' )
682          blx = INT( building_length_x / dx )
683          bly = INT( building_length_y / dy )
684          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
685          IF ( ABS( zw(bh  ) - building_height ) == &
686               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
687
688          IF ( building_wall_left == 9999999.9_wp )  THEN
689             building_wall_left = ( nx + 1 - blx ) / 2 * dx
690          ENDIF
691          bxl = INT ( building_wall_left / dx + 0.5_wp )
692          bxr = bxl + blx
693
694          IF ( building_wall_south == 9999999.9_wp )  THEN
695             building_wall_south = ( ny + 1 - bly ) / 2 * dy
696          ENDIF
697          bys = INT ( building_wall_south / dy + 0.5_wp )
698          byn = bys + bly
699
700          WRITE ( io, 271 )  building_length_x, building_length_y, &
701                             building_height, bxl, bxr, bys, byn
702
703       CASE ( 'single_street_canyon' )
704          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
705          IF ( ABS( zw(ch  ) - canyon_height ) == &
706               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
707          IF ( canyon_width_x /= 9999999.9_wp )  THEN
708!
709!--          Street canyon in y direction
710             cwx = NINT( canyon_width_x / dx )
711             IF ( canyon_wall_left == 9999999.9_wp )  THEN
712                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
713             ENDIF
714             cxl = NINT( canyon_wall_left / dx )
715             cxr = cxl + cwx
716             WRITE ( io, 272 )  'y', canyon_height, ch, 'u', cxl, cxr
717
718          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
719!
720!--          Street canyon in x direction
721             cwy = NINT( canyon_width_y / dy )
722             IF ( canyon_wall_south == 9999999.9_wp )  THEN
723                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
724             ENDIF
725             cys = NINT( canyon_wall_south / dy )
726             cyn = cys + cwy
727             WRITE ( io, 272 )  'x', canyon_height, ch, 'v', cys, cyn
728          ENDIF
729
730       CASE ( 'tunnel' )
731          IF ( tunnel_width_x /= 9999999.9_wp )  THEN
732!
733!--          Tunnel axis in y direction
734             IF ( tunnel_length == 9999999.9_wp  .OR.                          &
735                  tunnel_length >= ( nx + 1 ) * dx )  THEN
736                WRITE ( io, 273 )  'y', tunnel_height, tunnel_wall_depth,      &
737                                        tunnel_width_x
738             ELSE
739                WRITE ( io, 274 )  'y', tunnel_height, tunnel_wall_depth,      &
740                                        tunnel_width_x, tunnel_length
741             ENDIF
742
743          ELSEIF ( tunnel_width_y /= 9999999.9_wp )  THEN
744!
745!--          Tunnel axis in x direction
746             IF ( tunnel_length == 9999999.9_wp  .OR.                          &
747                  tunnel_length >= ( ny + 1 ) * dy )  THEN
748                WRITE ( io, 273 )  'x', tunnel_height, tunnel_wall_depth,      &
749                                        tunnel_width_y
750             ELSE
751                WRITE ( io, 274 )  'x', tunnel_height, tunnel_wall_depth,      &
752                                        tunnel_width_y, tunnel_length
753             ENDIF
754          ENDIF
755
756    END SELECT
757
758    IF ( TRIM( topography ) /= 'flat' )  THEN
759       IF ( TRIM( topography_grid_convention ) == ' ' )  THEN
760          IF ( TRIM( topography ) == 'single_building' .OR.  &
761               TRIM( topography ) == 'single_street_canyon' )  THEN
762             WRITE ( io, 278 )
763          ELSEIF ( TRIM( topography ) == 'read_from_file' )  THEN
764             WRITE ( io, 279 )
765          ENDIF
766       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_edge' )  THEN
767          WRITE ( io, 278 )
768       ELSEIF ( TRIM( topography_grid_convention ) == 'cell_center' )  THEN
769          WRITE ( io, 279 )
770       ENDIF
771    ENDIF
772
773!-- Complex terrain
774    IF ( complex_terrain )  THEN
775       WRITE( io, 280 ) 
776       IF ( turbulent_inflow )  THEN
777          WRITE( io, 281 )  zu(topo_top_ind(0,0,0))
778       ENDIF
779       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
780          WRITE( io, 282 )
781       ENDIF
782    ENDIF
783!
784!-- Boundary conditions
785    IF ( ibc_p_b == 0 )  THEN
786       r_lower = 'p(0)     = 0      |'
787    ELSEIF ( ibc_p_b == 1 )  THEN
788       r_lower = 'p(0)     = p(1)   |'
789    ENDIF
790    IF ( ibc_p_t == 0 )  THEN
791       r_upper  = 'p(nzt+1) = 0      |'
792    ELSE
793       r_upper  = 'p(nzt+1) = p(nzt) |'
794    ENDIF
795
796    IF ( ibc_uv_b == 0 )  THEN
797       r_lower = TRIM( r_lower ) // ' uv(0)     = -uv(1)                |'
798    ELSE
799       r_lower = TRIM( r_lower ) // ' uv(0)     = uv(1)                 |'
800    ENDIF
801    IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
802       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = 0                     |'
803    ELSEIF ( ibc_uv_t == 0 )  THEN
804       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = ug(nzt+1), vg(nzt+1)  |'
805    ELSE
806       r_upper  = TRIM( r_upper  ) // ' uv(nzt+1) = uv(nzt)               |'
807    ENDIF
808
809    IF ( ibc_pt_b == 0 )  THEN
810       IF ( land_surface )  THEN
811          r_lower = TRIM( r_lower ) // ' pt(0)     = from soil model'
812       ELSE
813          r_lower = TRIM( r_lower ) // ' pt(0)     = pt_surface'
814       ENDIF
815    ELSEIF ( ibc_pt_b == 1 )  THEN
816       r_lower = TRIM( r_lower ) // ' pt(0)     = pt(1)'
817    ELSEIF ( ibc_pt_b == 2 )  THEN
818       r_lower = TRIM( r_lower ) // ' pt(0)     = from coupled model'
819    ENDIF
820    IF ( ibc_pt_t == 0 )  THEN
821       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt_top'
822    ELSEIF( ibc_pt_t == 1 )  THEN
823       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt)'
824    ELSEIF( ibc_pt_t == 2 )  THEN
825       r_upper  = TRIM( r_upper  ) // ' pt(nzt+1) = pt(nzt) + dpt/dz_ini'
826
827    ENDIF
828
829    WRITE ( io, 300 )  r_lower, r_upper
830
831    IF ( .NOT. constant_diffusion )  THEN
832       IF ( ibc_e_b == 1 )  THEN
833          r_lower = 'e(0)     = e(1)'
834       ELSE
835          r_lower = 'e(0)     = e(1) = (u*/0.1)**2'
836       ENDIF
837       r_upper = 'e(nzt+1) = e(nzt) = e(nzt-1)'
838
839       WRITE ( io, 301 )  'e', r_lower, r_upper       
840
841    ENDIF
842
843    IF ( ocean_mode )  THEN
844       r_lower = 'sa(0)    = sa(1)'
845       IF ( ibc_sa_t == 0 )  THEN
846          r_upper =  'sa(nzt+1) = sa_surface'
847       ELSE
848          r_upper =  'sa(nzt+1) = sa(nzt)'
849       ENDIF
850       WRITE ( io, 301 ) 'sa', r_lower, r_upper
851    ENDIF
852
853    IF ( humidity )  THEN
854       IF ( ibc_q_b == 0 )  THEN
855          IF ( land_surface )  THEN
856             r_lower = 'q(0)     = from soil model'
857          ELSE
858             r_lower = 'q(0)     = q_surface'
859          ENDIF
860
861       ELSE
862          r_lower = 'q(0)      = q(1)'
863       ENDIF
864       IF ( ibc_q_t == 0 )  THEN
865          r_upper =  'q(nzt+1) = q_top'
866       ELSE
867          r_upper =  'q(nzt+1) = q(nzt) + dq/dz'
868       ENDIF
869       WRITE ( io, 301 ) 'q', r_lower, r_upper
870    ENDIF
871
872    IF ( passive_scalar )  THEN
873       IF ( ibc_s_b == 0 )  THEN
874          r_lower = 's(0)      = s_surface'
875       ELSE
876          r_lower = 's(0)      = s(1)'
877       ENDIF
878       IF ( ibc_s_t == 0 )  THEN
879          r_upper =  's(nzt+1) = s_top'
880       ELSEIF ( ibc_s_t == 1 )  THEN
881          r_upper =  's(nzt+1) = s(nzt)'
882       ELSEIF ( ibc_s_t == 2 )  THEN
883          r_upper =  's(nzt+1) = s(nzt) + ds/dz'
884       ENDIF
885       WRITE ( io, 301 ) 's', r_lower, r_upper
886    ENDIF
887
888    IF ( use_surface_fluxes )  THEN
889       WRITE ( io, 303 )
890       IF ( constant_heatflux )  THEN
891          IF ( large_scale_forcing .AND. lsf_surf )  THEN
892             IF ( surf_def_h(0)%ns >= 1 )  WRITE ( io, 306 )  surf_def_h(0)%shf(1)
893          ELSE
894             WRITE ( io, 306 )  surface_heatflux
895          ENDIF
896          IF ( random_heatflux )  WRITE ( io, 307 )
897       ENDIF
898       IF ( humidity  .AND.  constant_waterflux )  THEN
899          IF ( large_scale_forcing .AND. lsf_surf )  THEN
900             WRITE ( io, 311 ) surf_def_h(0)%qsws(1)
901          ELSE
902             WRITE ( io, 311 ) surface_waterflux
903          ENDIF
904       ENDIF
905       IF ( passive_scalar  .AND.  constant_scalarflux )  THEN
906          WRITE ( io, 313 ) surface_scalarflux
907       ENDIF
908    ENDIF
909
910    IF ( use_top_fluxes )  THEN
911       WRITE ( io, 304 )
912       IF ( coupling_mode == 'uncoupled' )  THEN
913          WRITE ( io, 320 )  top_momentumflux_u, top_momentumflux_v
914          IF ( constant_top_heatflux )  THEN
915             WRITE ( io, 306 )  top_heatflux
916          ENDIF
917       ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
918          WRITE ( io, 316 )
919       ENDIF
920       IF ( ocean_mode  .AND.  constant_top_salinityflux )                          &
921          WRITE ( io, 309 )  top_salinityflux
922       IF ( humidity       )  WRITE ( io, 315 )
923       IF ( passive_scalar .AND.  constant_top_scalarflux )                    &
924          WRITE ( io, 302 ) top_scalarflux
925    ENDIF
926
927    IF ( constant_flux_layer )  THEN
928       WRITE ( io, 305 )  (zu(1)-zu(0)), roughness_length,                     &
929                          z0h_factor*roughness_length, kappa,                  &
930                          zeta_min, zeta_max
931       IF ( .NOT. constant_heatflux )  WRITE ( io, 308 )
932       IF ( humidity  .AND.  .NOT. constant_waterflux )  THEN
933          WRITE ( io, 312 )
934       ENDIF
935       IF ( passive_scalar  .AND.  .NOT. constant_scalarflux )  THEN
936          WRITE ( io, 314 )
937       ENDIF
938    ELSE
939       IF ( INDEX(initializing_actions, 'set_1d-model_profiles') /= 0 )  THEN
940          WRITE ( io, 310 )  zeta_min, zeta_max
941       ENDIF
942    ENDIF
943
944    WRITE ( io, 317 )  bc_lr, bc_ns
945    IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
946       WRITE ( io, 318 )  use_cmax, pt_damping_width, pt_damping_factor       
947       IF ( turbulent_inflow )  THEN
948          IF ( y_shift == 0 ) THEN
949             WRITE ( io, 319 )  recycling_width, recycling_plane, &
950                                inflow_damping_height, inflow_damping_width
951          ELSE
952             WRITE ( io, 322 )  y_shift, recycling_width, recycling_plane, &
953                                inflow_damping_height, inflow_damping_width
954          END IF
955       ENDIF
956       IF ( turbulent_outflow )  THEN
957          WRITE ( io, 323 )  outflow_source_plane, INT(outflow_source_plane/dx)
958       ENDIF
959    ENDIF
960
961!
962!-- Initial Profiles
963    WRITE ( io, 321 )
964!
965!-- Initial wind profiles
966    IF ( u_profile(1) /= 9999999.9_wp )  WRITE ( io, 427 )
967
968!
969!-- Initial temperature profile
970!-- Building output strings, starting with surface temperature
971    WRITE ( temperatures, '(F6.2)' )  pt_surface
972    gradients = '------'
973    slices = '     0'
974    coordinates = '   0.0'
975    i = 1
976    DO  WHILE ( pt_vertical_gradient_level_ind(i) /= -9999 )
977
978       WRITE (coor_chr,'(F7.2)')  pt_init(pt_vertical_gradient_level_ind(i))
979       temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
980
981       WRITE (coor_chr,'(F7.2)')  pt_vertical_gradient(i)
982       gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
983
984       WRITE (coor_chr,'(I7)')  pt_vertical_gradient_level_ind(i)
985       slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
986
987       WRITE (coor_chr,'(F7.1)')  pt_vertical_gradient_level(i)
988       coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
989
990       IF ( i == 10 )  THEN
991          EXIT
992       ELSE
993          i = i + 1
994       ENDIF
995
996    ENDDO
997
998    IF ( .NOT. nudging )  THEN
999       WRITE ( io, 420 )  TRIM( coordinates ), TRIM( temperatures ), &
1000                          TRIM( gradients ), TRIM( slices )
1001    ELSE
1002       WRITE ( io, 428 ) 
1003    ENDIF
1004
1005!
1006!-- Initial humidity profile
1007!-- Building output strings, starting with surface humidity
1008    IF ( humidity )  THEN
1009       WRITE ( temperatures, '(E8.1)' )  q_surface
1010       gradients = '--------'
1011       slices = '       0'
1012       coordinates = '     0.0'
1013       i = 1
1014       DO  WHILE ( q_vertical_gradient_level_ind(i) /= -9999 )
1015         
1016          WRITE (coor_chr,'(E8.1,4X)')  q_init(q_vertical_gradient_level_ind(i))
1017          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1018
1019          WRITE (coor_chr,'(E8.1,4X)')  q_vertical_gradient(i)
1020          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1021         
1022          WRITE (coor_chr,'(I8,4X)')  q_vertical_gradient_level_ind(i)
1023          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1024         
1025          WRITE (coor_chr,'(F8.1,4X)')  q_vertical_gradient_level(i)
1026          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1027
1028          IF ( i == 10 )  THEN
1029             EXIT
1030          ELSE
1031             i = i + 1
1032          ENDIF
1033
1034       ENDDO
1035
1036       IF ( .NOT. nudging )  THEN
1037          WRITE ( io, 421 )  TRIM( coordinates ), TRIM( temperatures ),        &
1038                             TRIM( gradients ), TRIM( slices )
1039       ENDIF
1040    ENDIF
1041!
1042!-- Initial scalar profile
1043!-- Building output strings, starting with surface humidity
1044    IF ( passive_scalar )  THEN
1045       WRITE ( temperatures, '(E8.1)' )  s_surface
1046       gradients = '--------'
1047       slices = '       0'
1048       coordinates = '     0.0'
1049       i = 1
1050       DO  WHILE ( s_vertical_gradient_level_ind(i) /= -9999 )
1051         
1052          WRITE (coor_chr,'(E8.1,4X)')  s_init(s_vertical_gradient_level_ind(i))
1053          temperatures = TRIM( temperatures ) // '  ' // TRIM( coor_chr )
1054
1055          WRITE (coor_chr,'(E8.1,4X)')  s_vertical_gradient(i)
1056          gradients = TRIM( gradients ) // '  ' // TRIM( coor_chr )
1057         
1058          WRITE (coor_chr,'(I8,4X)')  s_vertical_gradient_level_ind(i)
1059          slices = TRIM( slices ) // '  ' // TRIM( coor_chr )
1060         
1061          WRITE (coor_chr,'(F8.1,4X)')  s_vertical_gradient_level(i)
1062          coordinates = TRIM( coordinates ) // '  '  // TRIM( coor_chr )
1063
1064          IF ( i == 10 )  THEN
1065             EXIT
1066          ELSE
1067             i = i + 1
1068          ENDIF
1069
1070       ENDDO
1071
1072       WRITE ( io, 422 )  TRIM( coordinates ), TRIM( temperatures ),           &
1073                          TRIM( gradients ), TRIM( slices )
1074    ENDIF   
1075
1076!
1077!-- Initial salinity profile
1078!-- Building output strings, starting with surface salinity
1079    IF ( ocean_mode )  THEN
1080       WRITE ( temperatures, '(F6.2)' )  sa_surface
1081       gradients = '------'
1082       slices = '     0'
1083       coordinates = '   0.0'
1084       i = 1
1085       DO  WHILE ( sa_vertical_gradient_level_ind(i) /= -9999 )
1086
1087          WRITE (coor_chr,'(F7.2)')  sa_init(sa_vertical_gradient_level_ind(i))
1088          temperatures = TRIM( temperatures ) // ' ' // TRIM( coor_chr )
1089
1090          WRITE (coor_chr,'(F7.2)')  sa_vertical_gradient(i)
1091          gradients = TRIM( gradients ) // ' ' // TRIM( coor_chr )
1092
1093          WRITE (coor_chr,'(I7)')  sa_vertical_gradient_level_ind(i)
1094          slices = TRIM( slices ) // ' ' // TRIM( coor_chr )
1095
1096          WRITE (coor_chr,'(F7.1)')  sa_vertical_gradient_level(i)
1097          coordinates = TRIM( coordinates ) // ' '  // TRIM( coor_chr )
1098
1099          IF ( i == 10 )  THEN
1100             EXIT
1101          ELSE
1102             i = i + 1
1103          ENDIF
1104
1105       ENDDO
1106
1107       WRITE ( io, 425 )  TRIM( coordinates ), TRIM( temperatures ), &
1108                          TRIM( gradients ), TRIM( slices )
1109    ENDIF
1110
1111
1112!
1113!-- Listing of 1D-profiles
1114    WRITE ( io, 325 )  dt_dopr_listing
1115    IF ( averaging_interval_pr /= 0.0_wp )  THEN
1116       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
1117    ENDIF
1118
1119!
1120!-- DATA output
1121    WRITE ( io, 330 )
1122    IF ( averaging_interval_pr /= 0.0_wp )  THEN
1123       WRITE ( io, 326 )  averaging_interval_pr, dt_averaging_input_pr
1124    ENDIF
1125
1126!
1127!-- 1D-profiles
1128    dopr_chr = 'Profile:'
1129    IF ( dopr_n /= 0 )  THEN
1130       WRITE ( io, 331 )
1131
1132       output_format = ''
1133       output_format = netcdf_data_format_string
1134       IF ( netcdf_deflate == 0 )  THEN
1135          WRITE ( io, 344 )  output_format
1136       ELSE
1137          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1138       ENDIF
1139
1140       DO  i = 1, dopr_n
1141          dopr_chr = TRIM( dopr_chr ) // ' ' // TRIM( data_output_pr(i) ) // ','
1142          IF ( LEN_TRIM( dopr_chr ) >= 60 )  THEN
1143             WRITE ( io, 332 )  dopr_chr
1144             dopr_chr = '       :'
1145          ENDIF
1146       ENDDO
1147
1148       IF ( dopr_chr /= '' )  THEN
1149          WRITE ( io, 332 )  dopr_chr
1150       ENDIF
1151       WRITE ( io, 333 )  dt_dopr, averaging_interval_pr, dt_averaging_input_pr
1152       IF ( skip_time_dopr /= 0.0_wp )  WRITE ( io, 339 )  skip_time_dopr
1153    ENDIF
1154
1155!
1156!-- 2D-arrays
1157    DO  av = 0, 1
1158
1159       i = 1
1160       do2d_xy = ''
1161       do2d_xz = ''
1162       do2d_yz = ''
1163       DO  WHILE ( do2d(av,i) /= ' ' )
1164
1165          l = MAX( 2, LEN_TRIM( do2d(av,i) ) )
1166          do2d_mode = do2d(av,i)(l-1:l)
1167
1168          SELECT CASE ( do2d_mode )
1169             CASE ( 'xy' )
1170                ll = LEN_TRIM( do2d_xy )
1171                do2d_xy = do2d_xy(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1172             CASE ( 'xz' )
1173                ll = LEN_TRIM( do2d_xz )
1174                do2d_xz = do2d_xz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1175             CASE ( 'yz' )
1176                ll = LEN_TRIM( do2d_yz )
1177                do2d_yz = do2d_yz(1:ll) // ' ' // do2d(av,i)(1:l-3) // ','
1178          END SELECT
1179
1180          i = i + 1
1181
1182       ENDDO
1183
1184       IF ( ( ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  .OR.    &
1185              ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  .OR.    &
1186              ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 ) ) )  THEN
1187
1188          IF (  av == 0 )  THEN
1189             WRITE ( io, 334 )  ''
1190          ELSE
1191             WRITE ( io, 334 )  '(time-averaged)'
1192          ENDIF
1193
1194          IF ( do2d_at_begin )  THEN
1195             begin_chr = 'and at the start'
1196          ELSE
1197             begin_chr = ''
1198          ENDIF
1199
1200          output_format = ''
1201          output_format = netcdf_data_format_string
1202          IF ( netcdf_deflate == 0 )  THEN
1203             WRITE ( io, 344 )  output_format
1204          ELSE
1205             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1206          ENDIF
1207
1208          IF ( do2d_xy /= ''  .AND.  section(1,1) /= -9999 )  THEN
1209             i = 1
1210             slices = '/'
1211             coordinates = '/'
1212!
1213!--          Building strings with index and coordinate information of the
1214!--          slices
1215             DO  WHILE ( section(i,1) /= -9999 )
1216
1217                WRITE (section_chr,'(I5)')  section(i,1)
1218                section_chr = ADJUSTL( section_chr )
1219                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1220
1221                IF ( section(i,1) == -1 )  THEN
1222                   WRITE (coor_chr,'(F10.1)')  -1.0_wp
1223                ELSE
1224                   WRITE (coor_chr,'(F10.1)')  zu(section(i,1))
1225                ENDIF
1226                coor_chr = ADJUSTL( coor_chr )
1227                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1228
1229                i = i + 1
1230             ENDDO
1231             IF ( av == 0 )  THEN
1232                WRITE ( io, 335 )  'XY', do2d_xy, dt_do2d_xy, &
1233                                   TRIM( begin_chr ), 'k', TRIM( slices ), &
1234                                   TRIM( coordinates )
1235                IF ( skip_time_do2d_xy /= 0.0_wp )  THEN
1236                   WRITE ( io, 339 )  skip_time_do2d_xy
1237                ENDIF
1238             ELSE
1239                WRITE ( io, 342 )  'XY', do2d_xy, dt_data_output_av, &
1240                                   TRIM( begin_chr ), averaging_interval, &
1241                                   dt_averaging_input, 'k', TRIM( slices ), &
1242                                   TRIM( coordinates )
1243                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
1244                   WRITE ( io, 339 )  skip_time_data_output_av
1245                ENDIF
1246             ENDIF
1247             IF ( netcdf_data_format > 4 )  THEN
1248                WRITE ( io, 352 )  ntdim_2d_xy(av)
1249             ELSE
1250                WRITE ( io, 353 )
1251             ENDIF
1252          ENDIF
1253
1254          IF ( do2d_xz /= ''  .AND.  section(1,2) /= -9999 )  THEN
1255             i = 1
1256             slices = '/'
1257             coordinates = '/'
1258!
1259!--          Building strings with index and coordinate information of the
1260!--          slices
1261             DO  WHILE ( section(i,2) /= -9999 )
1262
1263                WRITE (section_chr,'(I5)')  section(i,2)
1264                section_chr = ADJUSTL( section_chr )
1265                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1266
1267                WRITE (coor_chr,'(F10.1)')  section(i,2) * dy
1268                coor_chr = ADJUSTL( coor_chr )
1269                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1270
1271                i = i + 1
1272             ENDDO
1273             IF ( av == 0 )  THEN
1274                WRITE ( io, 335 )  'XZ', do2d_xz, dt_do2d_xz, &
1275                                   TRIM( begin_chr ), 'j', TRIM( slices ), &
1276                                   TRIM( coordinates )
1277                IF ( skip_time_do2d_xz /= 0.0_wp )  THEN
1278                   WRITE ( io, 339 )  skip_time_do2d_xz
1279                ENDIF
1280             ELSE
1281                WRITE ( io, 342 )  'XZ', do2d_xz, dt_data_output_av, &
1282                                   TRIM( begin_chr ), averaging_interval, &
1283                                   dt_averaging_input, 'j', TRIM( slices ), &
1284                                   TRIM( coordinates )
1285                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
1286                   WRITE ( io, 339 )  skip_time_data_output_av
1287                ENDIF
1288             ENDIF
1289             IF ( netcdf_data_format > 4 )  THEN
1290                WRITE ( io, 352 )  ntdim_2d_xz(av)
1291             ELSE
1292                WRITE ( io, 353 )
1293             ENDIF
1294          ENDIF
1295
1296          IF ( do2d_yz /= ''  .AND.  section(1,3) /= -9999 )  THEN
1297             i = 1
1298             slices = '/'
1299             coordinates = '/'
1300!
1301!--          Building strings with index and coordinate information of the
1302!--          slices
1303             DO  WHILE ( section(i,3) /= -9999 )
1304
1305                WRITE (section_chr,'(I5)')  section(i,3)
1306                section_chr = ADJUSTL( section_chr )
1307                slices = TRIM( slices ) // TRIM( section_chr ) // '/'
1308
1309                WRITE (coor_chr,'(F10.1)')  section(i,3) * dx
1310                coor_chr = ADJUSTL( coor_chr )
1311                coordinates = TRIM( coordinates ) // TRIM( coor_chr ) // '/'
1312
1313                i = i + 1
1314             ENDDO
1315             IF ( av == 0 )  THEN
1316                WRITE ( io, 335 )  'YZ', do2d_yz, dt_do2d_yz, &
1317                                   TRIM( begin_chr ), 'i', TRIM( slices ), &
1318                                   TRIM( coordinates )
1319                IF ( skip_time_do2d_yz /= 0.0_wp )  THEN
1320                   WRITE ( io, 339 )  skip_time_do2d_yz
1321                ENDIF
1322             ELSE
1323                WRITE ( io, 342 )  'YZ', do2d_yz, dt_data_output_av, &
1324                                   TRIM( begin_chr ), averaging_interval, &
1325                                   dt_averaging_input, 'i', TRIM( slices ), &
1326                                   TRIM( coordinates )
1327                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
1328                   WRITE ( io, 339 )  skip_time_data_output_av
1329                ENDIF
1330             ENDIF
1331             IF ( netcdf_data_format > 4 )  THEN
1332                WRITE ( io, 352 )  ntdim_2d_yz(av)
1333             ELSE
1334                WRITE ( io, 353 )
1335             ENDIF
1336          ENDIF
1337
1338       ENDIF
1339
1340    ENDDO
1341
1342!
1343!-- 3d-arrays
1344    DO  av = 0, 1
1345
1346       i = 1
1347       do3d_chr = ''
1348       DO  WHILE ( do3d(av,i) /= ' ' )
1349
1350          do3d_chr = TRIM( do3d_chr ) // ' ' // TRIM( do3d(av,i) ) // ','
1351          i = i + 1
1352
1353       ENDDO
1354
1355       IF ( do3d_chr /= '' )  THEN
1356          IF ( av == 0 )  THEN
1357             WRITE ( io, 336 )  ''
1358          ELSE
1359             WRITE ( io, 336 )  '(time-averaged)'
1360          ENDIF
1361
1362          output_format = netcdf_data_format_string
1363          IF ( netcdf_deflate == 0 )  THEN
1364             WRITE ( io, 344 )  output_format
1365          ELSE
1366             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1367          ENDIF
1368
1369          IF ( do3d_at_begin )  THEN
1370             begin_chr = 'and at the start'
1371          ELSE
1372             begin_chr = ''
1373          ENDIF
1374          IF ( av == 0 )  THEN
1375             WRITE ( io, 337 )  do3d_chr, dt_do3d, TRIM( begin_chr ), &
1376                                zu(nz_do3d), nz_do3d
1377          ELSE
1378             WRITE ( io, 343 )  do3d_chr, dt_data_output_av,           &
1379                                TRIM( begin_chr ), averaging_interval, &
1380                                dt_averaging_input, zu(nz_do3d), nz_do3d
1381          ENDIF
1382
1383          IF ( netcdf_data_format > 4 )  THEN
1384             WRITE ( io, 352 )  ntdim_3d(av)
1385          ELSE
1386             WRITE ( io, 353 )
1387          ENDIF
1388
1389          IF ( av == 0 )  THEN
1390             IF ( skip_time_do3d /= 0.0_wp )  THEN
1391                WRITE ( io, 339 )  skip_time_do3d
1392             ENDIF
1393          ELSE
1394             IF ( skip_time_data_output_av /= 0.0_wp )  THEN
1395                WRITE ( io, 339 )  skip_time_data_output_av
1396             ENDIF
1397          ENDIF
1398
1399       ENDIF
1400
1401    ENDDO
1402
1403!
1404!-- masked arrays
1405    IF ( masks > 0 )  WRITE ( io, 345 )  &
1406         mask_scale_x, mask_scale_y, mask_scale_z
1407    DO  mid = 1, masks
1408       DO  av = 0, 1
1409
1410          i = 1
1411          domask_chr = ''
1412          DO  WHILE ( domask(mid,av,i) /= ' ' )
1413             domask_chr = TRIM( domask_chr ) // ' ' //  &
1414                          TRIM( domask(mid,av,i) ) // ','
1415             i = i + 1
1416          ENDDO
1417
1418          IF ( domask_chr /= '' )  THEN
1419             IF ( av == 0 )  THEN
1420                WRITE ( io, 346 )  '', mid
1421             ELSE
1422                WRITE ( io, 346 )  ' (time-averaged)', mid
1423             ENDIF
1424
1425             output_format = netcdf_data_format_string
1426!--          Parallel output not implemented for mask data, hence
1427!--          output_format must be adjusted.
1428             IF ( netcdf_data_format == 5 ) output_format = 'netCDF4/HDF5'
1429             IF ( netcdf_data_format == 6 ) output_format = 'netCDF4/HDF5 classic'
1430             IF ( netcdf_deflate == 0 )  THEN
1431                WRITE ( io, 344 )  output_format
1432             ELSE
1433                WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1434             ENDIF
1435
1436             IF ( av == 0 )  THEN
1437                WRITE ( io, 347 )  domask_chr, dt_domask(mid)
1438             ELSE
1439                WRITE ( io, 348 )  domask_chr, dt_data_output_av, &
1440                                   averaging_interval, dt_averaging_input
1441             ENDIF
1442
1443             IF ( av == 0 )  THEN
1444                IF ( skip_time_domask(mid) /= 0.0_wp )  THEN
1445                   WRITE ( io, 339 )  skip_time_domask(mid)
1446                ENDIF
1447             ELSE
1448                IF ( skip_time_data_output_av /= 0.0_wp )  THEN
1449                   WRITE ( io, 339 )  skip_time_data_output_av
1450                ENDIF
1451             ENDIF
1452!
1453!--          output locations
1454             DO  dim = 1, 3
1455                IF ( mask(mid,dim,1) >= 0.0_wp )  THEN
1456                   count = 0
1457                   DO  WHILE ( mask(mid,dim,count+1) >= 0.0_wp )
1458                      count = count + 1
1459                   ENDDO
1460                   WRITE ( io, 349 )  dir(dim), dir(dim), mid, dir(dim), &
1461                                      mask(mid,dim,:count)
1462                ELSEIF ( mask_loop(mid,dim,1) < 0.0_wp .AND.  &
1463                         mask_loop(mid,dim,2) < 0.0_wp .AND.  &
1464                         mask_loop(mid,dim,3) == 0.0_wp )  THEN
1465                   WRITE ( io, 350 )  dir(dim), dir(dim)
1466                ELSEIF ( mask_loop(mid,dim,3) == 0.0_wp )  THEN
1467                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1468                                      mask_loop(mid,dim,1:2)
1469                ELSE
1470                   WRITE ( io, 351 )  dir(dim), dir(dim), mid, dir(dim), &
1471                                      mask_loop(mid,dim,1:3)
1472                ENDIF
1473             ENDDO
1474          ENDIF
1475
1476       ENDDO
1477    ENDDO
1478
1479!
1480!-- Timeseries
1481    IF ( dt_dots /= 9999999.9_wp )  THEN
1482       WRITE ( io, 340 )
1483
1484       output_format = netcdf_data_format_string
1485       IF ( netcdf_deflate == 0 )  THEN
1486          WRITE ( io, 344 )  output_format
1487       ELSE
1488          WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
1489       ENDIF
1490       WRITE ( io, 341 )  dt_dots
1491    ENDIF
1492
1493    WRITE ( io, 99 )
1494
1495!
1496!-- Physical quantities
1497    WRITE ( io, 400 )
1498
1499!
1500!-- Geostrophic parameters
1501    WRITE ( io, 410 )  latitude, longitude, rotation_angle, omega, f, fs
1502
1503!
1504!-- Day and time during model start
1505    CALL get_date_time( 0.0_wp, date_time_str=date_time_str )
1506    WRITE ( io, 456 )  TRIM( date_time_str )
1507
1508!
1509!-- Other quantities
1510    WRITE ( io, 411 )  g
1511
1512    WRITE ( io, 412 )  TRIM( reference_state )
1513    IF ( use_single_reference_value )  THEN
1514       IF ( ocean_mode )  THEN
1515          WRITE ( io, 413 )  prho_reference
1516       ELSE
1517          WRITE ( io, 414 )  pt_reference
1518       ENDIF
1519    ENDIF
1520
1521!
1522!-- Cloud physcis parameters / quantities / numerical methods
1523    WRITE ( io, 430 )
1524    IF ( humidity .AND. .NOT. bulk_cloud_model .AND. .NOT. cloud_droplets)  THEN
1525       WRITE ( io, 431 )
1526    ENDIF
1527!
1528!-- LES / turbulence parameters
1529    WRITE ( io, 450 )
1530
1531!--
1532! ... LES-constants used must still be added here
1533!--
1534    IF ( constant_diffusion )  THEN
1535       WRITE ( io, 451 )  km_constant, km_constant/prandtl_number, &
1536                          prandtl_number
1537    ENDIF
1538    IF ( .NOT. constant_diffusion)  THEN
1539       IF ( e_init > 0.0_wp )  WRITE ( io, 455 )  e_init
1540       IF ( e_min > 0.0_wp )  WRITE ( io, 454 )  e_min
1541       IF ( wall_adjustment )  WRITE ( io, 453 )  wall_adjustment_factor
1542    ENDIF
1543    IF ( rans_mode )  THEN
1544       WRITE ( io, 457 )  rans_const_c, rans_const_sigma
1545    ENDIF
1546!
1547!-- Special actions during the run
1548    WRITE ( io, 470 )
1549    IF ( create_disturbances )  THEN
1550       WRITE ( io, 471 )  dt_disturb, disturbance_amplitude,                   &
1551                          zu(disturbance_level_ind_b), disturbance_level_ind_b,&
1552                          zu(disturbance_level_ind_t), disturbance_level_ind_t
1553       IF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
1554          WRITE ( io, 472 )  inflow_disturbance_begin, inflow_disturbance_end
1555       ELSE
1556          WRITE ( io, 473 )  disturbance_energy_limit
1557       ENDIF
1558       WRITE ( io, 474 )  TRIM( random_generator )
1559    ENDIF
1560    IF ( pt_surface_initial_change /= 0.0_wp )  THEN
1561       WRITE ( io, 475 )  pt_surface_initial_change
1562    ENDIF
1563    IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
1564       WRITE ( io, 476 )  q_surface_initial_change       
1565    ENDIF
1566    IF ( passive_scalar  .AND.  q_surface_initial_change /= 0.0_wp )  THEN
1567       WRITE ( io, 477 )  q_surface_initial_change       
1568    ENDIF
1569
1570!
1571!-- Parameters of 1D-model
1572    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1573       WRITE ( io, 500 )  end_time_1d, dt_run_control_1d, dt_pr_1d, &
1574                          mixing_length_1d, dissipation_1d
1575       IF ( damp_level_ind_1d /= nzt+1 )  THEN
1576          WRITE ( io, 502 )  zu(damp_level_ind_1d), damp_level_ind_1d
1577       ENDIF
1578    ENDIF
1579
1580!
1581!-- Header information from other modules
1582    CALL module_interface_header( io )
1583
1584
1585    WRITE ( io, 99 )
1586
1587!
1588!-- Write buffer contents to disc immediately
1589    FLUSH( io )
1590
1591!
1592!-- Here the FORMATs start
1593
1594 99 FORMAT (1X,78('-'))
1595100 FORMAT (/1X,'******************************',4X,44('-')/        &
1596            1X,'* ',A,' *',4X,A/                               &
1597            1X,'******************************',4X,44('-'))
1598101 FORMAT (35X,'coupled run: ',A/ &
1599            35X,42('-'))
1600102 FORMAT (/' Date:               ',A10,4X,'Run:       ',A34/      &
1601            ' Time:                 ',A8,4X,'Run-No.:   ',I2.2/     &
1602            ' Run on host:        ',A10)
1603#if defined( __parallel )
1604103 FORMAT (' Number of PEs:',10X,I6,4X,'Processor grid (x,y): (',I4,',',I4, &
1605              ')',1X,A)
1606104 FORMAT (' Number of PEs:',10X,I6,4X,'Tasks:',I4,'   threads per task:',I4/ &
1607              35X,'Processor grid (x,y): (',I4,',',I4,')',1X,A)
1608107 FORMAT (35X,'A 1d-decomposition along ',A,' is used')
1609108 FORMAT (35X,'Max. # of parallel I/O streams is ',I5)
1610109 FORMAT (35X,'Precursor run for coupled atmos-ocean run'/ &
1611            35X,42('-'))
1612114 FORMAT (35X,'Coupled atmosphere-ocean run following'/ &
1613            35X,'independent precursor runs'/             &
1614            35X,42('-'))
1615#endif
1616110 FORMAT (/' Numerical Schemes:'/ &
1617             ' -----------------'/)
1618124 FORMAT (' --> Use the ',A,' turbulence closure (',A,' mode).')
1619121 FORMAT (' --> Use the ',A,' approximation for the model equations.')
1620111 FORMAT (' --> Solve perturbation pressure via FFT using ',A,' routines')
1621112 FORMAT (' --> Solve perturbation pressure via SOR-Red/Black-Schema'/ &
1622            '     Iterations (initial/other): ',I3,'/',I3,'  omega =',F6.3)
1623113 FORMAT (' --> Momentum advection via Piascek-Williams-Scheme (Form C3)', &
1624                  ' or Upstream')
1625115 FORMAT ('     FFT and transpositions are overlapping')
1626116 FORMAT (' --> Scalar advection via Piascek-Williams-Scheme (Form C3)', &
1627                  ' or Upstream')
1628118 FORMAT (' --> Scalar advection via Bott-Chlond-Scheme')
1629119 FORMAT (' --> Galilei-Transform applied to horizontal advection:'/ &
1630            '     translation velocity = ',A/ &
1631            '     distance advected ',A,':  ',F8.3,' km(x)  ',F8.3,' km(y)')
1632122 FORMAT (' --> Time differencing scheme: ',A)
1633123 FORMAT (' --> Rayleigh-Damping active, starts ',A,' z = ',F8.2,' m'/ &
1634            '     maximum damping coefficient:',F6.3, ' 1/s')
1635129 FORMAT (' --> Additional prognostic equation for the specific humidity')
1636130 FORMAT (' --> Additional prognostic equation for the total water content')
1637131 FORMAT (' --> No pt-equation solved. Neutral stratification with pt = ', &
1638                  F6.2, ' K assumed')
1639134 FORMAT (' --> Additional prognostic equation for a passive scalar')
1640135 FORMAT (' --> Solve perturbation pressure via ',A,' method (', &
1641                  A,'-cycle)'/ &
1642            '     number of grid levels:                   ',I2/ &
1643            '     Gauss-Seidel red/black iterations:       ',I2)
1644136 FORMAT ('     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1645                  I3,')')
1646137 FORMAT ('     level data gathered on PE0 at level:     ',I2/ &
1647            '     gridpoints of coarsest subdomain (x,y,z): (',I3,',',I3,',', &
1648                  I3,')'/ &
1649            '     gridpoints of coarsest domain (x,y,z):    (',I3,',',I3,',', &
1650                  I3,')')
1651139 FORMAT (' --> Loop optimization method: ',A)
1652140 FORMAT ('     maximum residual allowed:                ',E10.3)
1653141 FORMAT ('     fixed number of multigrid cycles:        ',I4)
1654142 FORMAT ('     perturbation pressure is calculated at every Runge-Kutta ', &
1655                  'step')
1656143 FORMAT ('     Euler/upstream scheme is used for the SGS turbulent ', &
1657                  'kinetic energy')
1658144 FORMAT ('     masking method is used')
1659150 FORMAT (' --> Volume flow at the right and north boundary will be ', &
1660                  'conserved'/ &
1661            '     using the ',A,' mode')
1662151 FORMAT ('     with u_bulk = ',F7.3,' m/s and v_bulk = ',F7.3,' m/s')
1663152 FORMAT (' --> External pressure gradient directly prescribed by the user:',&
1664           /'     ',2(1X,E12.5),'Pa/m in x/y direction', &
1665           /'     starting from dp_level_b =', F8.3, 'm', A /)
1666200 FORMAT (//' Run time and time step information:'/ &
1667             ' ----------------------------------'/)
1668201 FORMAT ( ' Timestep:             variable     maximum value: ',F6.3,' s', &
1669             '    CFL-factor:',F5.2)
1670202 FORMAT ( ' Timestep:          dt = ',F6.3,' s'/)
1671203 FORMAT ( ' Start time:        ',F11.3,' s'/ &
1672             ' End time:          ',F11.3,' s')
1673204 FORMAT ( A,F11.3,' s')
1674205 FORMAT ( A,F11.3,' s',5X,'restart every',17X,F11.3,' s')
1675206 FORMAT (/' Time reached:      ',F11.3,' s'/ &
1676             ' CPU-time used:       ',F9.3,' s     per timestep:                 ',F9.3,' s'/ &
1677             '                                      per second of simulated time: ',F9.3,' s')
1678207 FORMAT ( ' Spinup time:       ',F11.3,' s')
1679250 FORMAT (//' Computational grid and domain size:'/ &
1680              ' ----------------------------------'// &
1681              ' Grid length:      dx =    ',F8.3,' m    dy =    ',F8.3, ' m')
1682251 FORMAT (  /' Domain size:       x = ',F10.3,' m     y = ',F10.3, &
1683              ' m  z(u) = ',F10.3,' m'/)
1684253 FORMAT ( '                dz(',I1,') =    ', F8.3, ' m')
1685254 FORMAT (//' Number of gridpoints (x,y,z):  (0:',I4,', 0:',I4,', 0:',I4,')'/ &
1686            ' Subdomain size (x,y,z):        (  ',I4,',   ',I4,',   ',I4,')'/)
1687260 FORMAT (/' The model has a slope in x-direction. Inclination angle: ',F6.2,&
1688             ' degrees')
1689270 FORMAT (//' Topography information:'/ &
1690              ' ----------------------'// &
1691              1X,'Topography: ',A)
1692271 FORMAT (  ' Building size (x/y/z) in m: ',F5.1,' / ',F5.1,' / ',F5.1/ &
1693              ' Horizontal index bounds (l/r/s/n): ',I4,' / ',I4,' / ',I4, &
1694                ' / ',I4)
1695272 FORMAT (  ' Single quasi-2D street canyon of infinite length in ',A, &
1696              ' direction' / &
1697              ' Canyon height: ', F6.2, 'm, ch = ', I4, '.'      / &
1698              ' Canyon position (',A,'-walls): cxl = ', I4,', cxr = ', I4, '.')
1699273 FORMAT (  ' Tunnel of infinite length in ',A, &
1700              ' direction' / &
1701              ' Tunnel height: ', F6.2, / &
1702              ' Tunnel-wall depth: ', F6.2      / &
1703              ' Tunnel width: ', F6.2 )
1704274 FORMAT (  ' Tunnel in ', A, ' direction.' / &
1705              ' Tunnel height: ', F6.2, / &   
1706              ' Tunnel-wall depth: ', F6.2      / &
1707              ' Tunnel width: ', F6.2, / &
1708              ' Tunnel length: ', F6.2 )
1709278 FORMAT (' Topography grid definition convention:'/ &
1710            ' cell edge (staggered grid points'/  &
1711            ' (u in x-direction, v in y-direction))' /)
1712279 FORMAT (' Topography grid definition convention:'/ &
1713            ' cell center (scalar grid points)' /)
1714280 FORMAT (' Complex terrain simulation is activated.')
1715281 FORMAT ('    --> Mean inflow profiles are adjusted.' / &
1716            '    --> Elevation of inflow boundary: ', F7.1, ' m' )
1717282 FORMAT ('    --> Initial data from 3D-precursor run is shifted' / &
1718            '        vertically depending on local surface height.')
1719300 FORMAT (//' Boundary conditions:'/ &
1720             ' -------------------'// &
1721             '                     p                    uv             ', &
1722             '                     pt'// &
1723             ' B. bound.: ',A/ &
1724             ' T. bound.: ',A)
1725301 FORMAT (/'                     ',A// &
1726             ' B. bound.: ',A/ &
1727             ' T. bound.: ',A)
1728303 FORMAT (/' Bottom surface fluxes are used in diffusion terms at k=1')
1729304 FORMAT (/' Top surface fluxes are used in diffusion terms at k=nzt')
1730305 FORMAT (//'    Constant flux layer between bottom surface and first ',     &
1731              'computational u,v-level:'// &
1732             '       z_mo = ',F6.2,' m   z0 =',F7.4,' m   z0h =',F8.5,&
1733             ' m   kappa =',F5.2/ &
1734             '       Rif value range:   ',F8.2,' <= rif <=',F6.2)
1735306 FORMAT ('       Predefined constant heatflux:   ',F9.6,' K m/s')
1736307 FORMAT ('       Heatflux has a random normal distribution')
1737308 FORMAT ('       Predefined surface temperature')
1738309 FORMAT ('       Predefined constant salinityflux:   ',F9.6,' psu m/s')
1739310 FORMAT (//'    1D-Model:'// &
1740             '       Rif value range:   ',F6.2,' <= rif <=',F6.2)
1741311 FORMAT ('       Predefined constant humidity flux: ',E10.3,' kg/kg m/s')
1742312 FORMAT ('       Predefined surface humidity')
1743313 FORMAT ('       Predefined constant scalar flux: ',E10.3,' kg/(m**2 s)')
1744314 FORMAT ('       Predefined scalar value at the surface')
1745302 FORMAT ('       Predefined constant scalarflux:   ',F9.6,' kg/(m**2 s)')
1746315 FORMAT ('       Humidity flux at top surface is 0.0')
1747316 FORMAT ('       Sensible heatflux and momentum flux from coupled ', &
1748                    'atmosphere model')
1749317 FORMAT (//' Lateral boundaries:'/ &
1750            '       left/right:  ',A/    &
1751            '       north/south: ',A)
1752318 FORMAT (/'       use_cmax: ',L1 / &
1753            '       pt damping layer width = ',F8.2,' m, pt ', &
1754                    'damping factor =',F7.4)
1755319 FORMAT ('       turbulence recycling at inflow switched on'/ &
1756            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1757            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m')
1758320 FORMAT ('       Predefined constant momentumflux:  u: ',F9.6,' m**2/s**2'/ &
1759            '                                          v: ',F9.6,' m**2/s**2')
1760321 FORMAT (//' Initial profiles:'/ &
1761              ' ----------------')
1762322 FORMAT ('       turbulence recycling at inflow switched on'/ &
1763            '       y-shift of the recycled inflow turbulence is',I3,' PE'/ &
1764            '       width of recycling domain: ',F7.1,' m   grid index: ',I4/ &
1765            '       inflow damping height: ',F6.1,' m   width: ',F6.1,' m'/)
1766323 FORMAT ('       turbulent outflow conditon switched on'/ &
1767            '       position of outflow source plane: ',F7.1,' m   ', &
1768                    'grid index: ', I4)
1769325 FORMAT (//' List output:'/ &
1770             ' -----------'//  &
1771            '    1D-Profiles:'/    &
1772            '       Output every             ',F10.2,' s')
1773326 FORMAT ('       Time averaged over       ',F8.2,' s'/ &
1774            '       Averaging input every    ',F8.2,' s')
1775330 FORMAT (//' Data output:'/ &
1776             ' -----------'/)
1777331 FORMAT (/'    1D-Profiles:')
1778332 FORMAT (/'       ',A)
1779333 FORMAT ('       Output every             ',F8.2,' s',/ &
1780            '       Time averaged over       ',F8.2,' s'/ &
1781            '       Averaging input every    ',F8.2,' s')
1782334 FORMAT (/'    2D-Arrays',A,':')
1783335 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1784            '       Output every             ',F8.2,' s  ',A/ &
1785            '       Cross sections at ',A1,' = ',A/ &
1786            '       scalar-coordinates:   ',A,' m'/)
1787336 FORMAT (/'    3D-Arrays',A,':')
1788337 FORMAT (/'       Arrays: ',A/ &
1789            '       Output every             ',F8.2,' s  ',A/ &
1790            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1791339 FORMAT ('       No output during initial ',F8.2,' s')
1792340 FORMAT (/'    Time series:')
1793341 FORMAT ('       Output every             ',F8.2,' s'/)
1794342 FORMAT (/'       ',A2,'-cross-section  Arrays: ',A/ &
1795            '       Output every             ',F8.2,' s  ',A/ &
1796            '       Time averaged over       ',F8.2,' s'/ &
1797            '       Averaging input every    ',F8.2,' s'/ &
1798            '       Cross sections at ',A1,' = ',A/ &
1799            '       scalar-coordinates:   ',A,' m'/)
1800343 FORMAT (/'       Arrays: ',A/ &
1801            '       Output every             ',F8.2,' s  ',A/ &
1802            '       Time averaged over       ',F8.2,' s'/ &
1803            '       Averaging input every    ',F8.2,' s'/ &
1804            '       Upper output limit at    ',F8.2,' m  (GP ',I4,')'/)
1805344 FORMAT ('       Output format: ',A/)
1806345 FORMAT (/'    Scaling lengths for output locations of all subsequent mask IDs:',/ &
1807            '       mask_scale_x (in x-direction): ',F9.3, ' m',/ &
1808            '       mask_scale_y (in y-direction): ',F9.3, ' m',/ &
1809            '       mask_scale_z (in z-direction): ',F9.3, ' m' )
1810346 FORMAT (/'    Masked data output',A,' for mask ID ',I2, ':')
1811347 FORMAT ('       Variables: ',A/ &
1812            '       Output every             ',F8.2,' s')
1813348 FORMAT ('       Variables: ',A/ &
1814            '       Output every             ',F8.2,' s'/ &
1815            '       Time averaged over       ',F8.2,' s'/ &
1816            '       Averaging input every    ',F8.2,' s')
1817349 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1818            'mask_scale_',A,' predefined by array mask_',I2.2,'_',A,':'/ &
1819            13('       ',8(F8.2,',')/) )
1820350 FORMAT (/'       Output locations in ',A,'-direction: ', &
1821            'all gridpoints along ',A,'-direction (default).' )
1822351 FORMAT (/'       Output locations in ',A,'-direction in multiples of ', &
1823            'mask_scale_',A,' constructed from array mask_',I2.2,'_',A,'_loop:'/ &
1824            '          loop begin:',F8.2,', end:',F8.2,', stride:',F8.2 )
1825352 FORMAT  (/'       Number of output time levels allowed: ',I3 /)
1826353 FORMAT  (/'       Number of output time levels allowed: unlimited' /)
1827354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
1828400 FORMAT (//' Physical quantities:'/ &
1829              ' -------------------'/)
1830410 FORMAT ('    Geograph. latitude  :   latitude  = ',F5.1,' degr'/   &
1831            '    Geograph. longitude :   longitude = ',F5.1,' degr'/   &
1832            '    Rotation angle      :   rotation_angle = ',F5.1,' degr'/   &
1833            '    Angular velocity    :   omega  =',E10.3,' rad/s'/  &
1834            '    Coriolis parameter  :   f      = ',F9.6,' 1/s'/    &
1835            '                            f*     = ',F9.6,' 1/s')
1836411 FORMAT (/'    Gravity             :   g      = ',F4.1,' m/s**2')
1837412 FORMAT (/'    Reference state used in buoyancy terms: ',A)
1838413 FORMAT ('       Reference density in buoyancy terms: ',F8.3,' kg/m**3')
1839414 FORMAT ('       Reference temperature in buoyancy terms: ',F8.4,' K')
1840420 FORMAT (/'    Characteristic levels of the initial temperature profile:'// &
1841            '       Height:        ',A,'  m'/ &
1842            '       Temperature:   ',A,'  K'/ &
1843            '       Gradient:      ',A,'  K/100m'/ &
1844            '       Gridpoint:     ',A)
1845421 FORMAT (/'    Characteristic levels of the initial humidity profile:'// &
1846            '       Height:      ',A,'  m'/ &
1847            '       Humidity:    ',A,'  kg/kg'/ &
1848            '       Gradient:    ',A,'  (kg/kg)/100m'/ &
1849            '       Gridpoint:   ',A)
1850422 FORMAT (/'    Characteristic levels of the initial scalar profile:'// &
1851            '       Height:                  ',A,'  m'/ &
1852            '       Scalar concentration:    ',A,'  kg/m**3'/ &
1853            '       Gradient:                ',A,'  (kg/m**3)/100m'/ &
1854            '       Gridpoint:               ',A)
1855423 FORMAT (/'    Characteristic levels of the geo. wind component ug:'// &
1856            '       Height:      ',A,'  m'/ &
1857            '       ug:          ',A,'  m/s'/ &
1858            '       Gradient:    ',A,'  1/100s'/ &
1859            '       Gridpoint:   ',A)
1860424 FORMAT (/'    Characteristic levels of the geo. wind component vg:'// &
1861            '       Height:      ',A,'  m'/ &
1862            '       vg:          ',A,'  m/s'/ &
1863            '       Gradient:    ',A,'  1/100s'/ &
1864            '       Gridpoint:   ',A)
1865425 FORMAT (/'    Characteristic levels of the initial salinity profile:'// &
1866            '       Height:     ',A,'  m'/ &
1867            '       Salinity:   ',A,'  psu'/ &
1868            '       Gradient:   ',A,'  psu/100m'/ &
1869            '       Gridpoint:  ',A)
1870426 FORMAT (/'    Characteristic levels of the subsidence/ascent profile:'// &
1871            '       Height:      ',A,'  m'/ &
1872            '       w_subs:      ',A,'  m/s'/ &
1873            '       Gradient:    ',A,'  (m/s)/100m'/ &
1874            '       Gridpoint:   ',A)
1875427 FORMAT (/'    Initial wind profiles (u,v) are interpolated from given'// &
1876                  ' profiles')
1877428 FORMAT (/'    Initial profiles (u, v, pt, q) are taken from file '/ &
1878             '    NUDGING_DATA')
1879430 FORMAT (//' Cloud physics quantities / methods:'/ &
1880              ' ----------------------------------'/)
1881431 FORMAT ('    Humidity is considered, bu no condensation')
1882450 FORMAT (//' LES / Turbulence quantities:'/ &
1883              ' ---------------------------'/)
1884451 FORMAT ('    Diffusion coefficients are constant:'/ &
1885            '    Km = ',F6.2,' m**2/s   Kh = ',F6.2,' m**2/s   Pr = ',F5.2)
1886453 FORMAT ('    Mixing length is limited to',F5.2,' * z')
1887454 FORMAT ('    TKE is not allowed to fall below ',E9.2,' (m/s)**2')
1888455 FORMAT ('    initial TKE is prescribed as ',E9.2,' (m/s)**2')
1889456 FORMAT (/'    Date and time at model start : ',A)
1890457 FORMAT ('    RANS-mode constants: c_0 = ',F9.5/         &
1891            '                         c_1 = ',F9.5/         &
1892            '                         c_2 = ',F9.5/         &
1893            '                         c_3 = ',F9.5/         &
1894            '                         c_4 = ',F9.5/         &
1895            '                         sigma_e    = ',F9.5/  &
1896            '                         sigma_diss = ',F9.5)
1897470 FORMAT (//' Actions during the simulation:'/ &
1898              ' -----------------------------'/)
1899471 FORMAT ('    Disturbance impulse (u,v) every :   ',F6.2,' s'/            &
1900            '    Disturbance amplitude           :    ',F5.2, ' m/s'/       &
1901            '    Lower disturbance level         : ',F8.2,' m (GP ',I4,')'/  &
1902            '    Upper disturbance level         : ',F8.2,' m (GP ',I4,')')
1903472 FORMAT ('    Disturbances continued during the run from i/j =',I4, &
1904                 ' to i/j =',I4)
1905473 FORMAT ('    Disturbances cease as soon as the disturbance energy exceeds',&
1906                 F6.3, ' m**2/s**2')
1907474 FORMAT ('    Random number generator used    : ',A/)
1908475 FORMAT ('    The surface temperature is increased (or decreased, ', &
1909                 'respectively, if'/ &
1910            '    the value is negative) by ',F5.2,' K at the beginning of the',&
1911                 ' 3D-simulation'/)
1912476 FORMAT ('    The surface humidity is increased (or decreased, ',&
1913                 'respectively, if the'/ &
1914            '    value is negative) by ',E8.1,' kg/kg at the beginning of', &
1915                 ' the 3D-simulation'/)
1916477 FORMAT ('    The scalar value is increased at the surface (or decreased, ',&
1917                 'respectively, if the'/ &
1918            '    value is negative) by ',E8.1,' kg/m**3 at the beginning of', &
1919                 ' the 3D-simulation'/)
1920500 FORMAT (//' 1D-Model parameters:'/                           &
1921              ' -------------------'//                           &
1922            '    Simulation time:                   ',F8.1,' s'/ &
1923            '    Run-controll output every:         ',F8.1,' s'/ &
1924            '    Vertical profile output every:     ',F8.1,' s'/ &
1925            '    Mixing length calculation:         ',A/         &
1926            '    Dissipation calculation:           ',A/)
1927502 FORMAT ('    Damping layer starts from ',F7.1,' m (GP ',I4,')'/)
1928503 FORMAT (' --> Momentum advection via Wicker-Skamarock-Scheme 5th order')
1929504 FORMAT (' --> Scalar advection via Wicker-Skamarock-Scheme 5th order')
1930512 FORMAT (/' Date:               ',A10,6X,'Run:       ',A34/      &
1931            ' Time:                 ',A8,6X,'Run-No.:   ',I2.2/     &
1932            ' Run on host:        ',A10,6X,'En-No.:    ',I2.2)
1933#if defined( __parallel )
1934600 FORMAT (/' Nesting informations:'/ &
1935            ' --------------------'/ &
1936            ' Nesting mode:                     ',A/ &
1937            ' Nesting-datatransfer mode:        ',A// &
1938            ' Nest id  parent  number   lower left coordinates   name'/ &
1939            ' (*=me)     id    of PEs      x (m)     y (m)' )
1940601 FORMAT (2X,A1,1X,I2.2,6X,I2.2,5X,I5,5X,F8.2,2X,F8.2,5X,A)
1941#endif
1942
1943 END SUBROUTINE header
Note: See TracBrowser for help on using the repository browser.