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

Last change on this file since 4356 was 4339, checked in by suehring, 5 years ago

Bugfix, character strenght too short, caused crash on NEC

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