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

Last change on this file since 4304 was 4301, checked in by oliver.maas, 5 years ago

Deleted parameter recycling_yshift. y-shift in case of non-cyclic boundary conditions and turbulent_inflow = .TRUE. is now steered by parameter y_shift, that is also used in case of cyclic boundary conditions.

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