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

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

Consider rotation of model domain for claculation of Coriolis force:

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