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

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