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

Last change on this file since 4239 was 4227, checked in by gronemeier, 5 years ago

implement new palm_date_time_mod; replaced namelist parameters time_utc_init and day_of_year_init by origin_date_time

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