source: palm/trunk/SOURCE/check_parameters.f90 @ 4472

Last change on this file since 4472 was 4472, checked in by Giersch, 4 years ago

Profile output of the Kolmogorov length scale added

  • Property svn:keywords set to Id
  • Property svn:mergeinfo set to False
    /palm/branches/chemistry/SOURCE/check_parameters.f902047-3190,​3218-3297
    /palm/branches/forwind/SOURCE/check_parameters.f901564-1913
    /palm/branches/mosaik_M2/check_parameters.f902360-3471
    /palm/branches/palm4u/SOURCE/check_parameters.f902540-2692
    /palm/branches/rans/SOURCE/check_parameters.f902078-3128
    /palm/branches/resler/SOURCE/check_parameters.f902023-4391
    /palm/branches/salsa/SOURCE/check_parameters.f902503-3581
File size: 141.1 KB
Line 
1!> @file check_parameters.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-2020 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: check_parameters.f90 4472 2020-03-24 12:21:00Z Giersch $
27! Kolmogorov length scale eta added to profile output
28!
29! 4444 2020-03-05 15:59:50Z raasch
30! bugfix: cpp-directives for serial mode added
31!
32! 4392 2020-01-31 16:14:57Z pavelkrc
33! Some error numbers revised to prevent double usage
34!
35! 11:55:33Z oliver.maas
36! Checks for closed channel flow implemented
37!
38! 11:55:33Z oliver.maas
39! Move 2-m potential temperature output to diagnostic_output_quantities
40!
41! 11:55:33Z oliver.maas
42! removed message PA0421, concerning old parameter recycling_yshift
43!
44! 11:55:33Z oliver.maas
45! adjust message to the modified parameter recycling_yshift
46!
47! 11:55:33Z oliver.maas
48! Check if a cross section is specified if any output cross-section quantity
49! is given
50!
51! 11:55:33Z oliver.maas
52! Overwrite rotation_angle from namelist by value from static driver
53!
54! 11:55:33Z oliver.maas
55! removed conversion from recycle_absolute_quantities to raq, added check and
56! error message for correct input of recycling_method_for_thermodynamic_quantities
57!
58! 11:55:33Z oliver.maas
59! Corrected "Former revisions" section
60!
61! 11:55:33Z oliver.maas
62! bugfix error message: replaced PA184 by PA0184
63!
64! 11:55:33Z oliver.maas
65! added conversion from recycle_absolute_quantities to raq for recycling of
66! absolute quantities and added error message PA184 for not implemented quantities
67!
68! 4142 2019-08-05 12:38:31Z suehring
69! Consider spinup in number of output timesteps for averaged 2D output (merge
70! from branch resler).
71!
72! 4069 2019-07-01 14:05:51Z Giersch
73! Masked output running index mid has been introduced as a local variable to
74! avoid runtime error (Loop variable has been modified) in time_integration
75!
76! 4048 2019-06-21 21:00:21Z knoop
77! Moved tcm_check_data_output to module_interface
78!
79! 4039 2019-06-18 10:32:41Z suehring
80! Modularize diagnostic output
81!
82! 4017 2019-06-06 12:16:46Z schwenkel
83! output of turbulence intensity added
84!
85! 3933 2019-04-25 12:33:20Z kanani
86! Alphabetical resorting in CASE, condense settings for theta_2m* into one IF clause
87!
88! 3885 2019-04-11 11:29:34Z kanani
89! Changes related to global restructuring of location messages and introduction
90! of additional debug messages
91!
92! 3766 2019-02-26 16:23:41Z raasch
93! trim added to avoid truncation compiler warnings
94!
95! 3761 2019-02-25 15:31:42Z raasch
96! unused variables removed
97!
98! 3735 2019-02-12 09:52:40Z dom_dwd_user
99! Passing variable j (averaged output?) to
100! module_interface.f90:chem_check_data_output.
101!
102! 3705 2019-01-29 19:56:39Z suehring
103! bugfix: renamed thetav_t to vtheta_t
104!
105! 3702 2019-01-28 13:19:30Z gronemeier
106! most_method removed
107!
108! 3655 2019-01-07 16:51:22Z knoop
109! Formatting
110!
111! Revision 1.1  1997/08/26 06:29:23  raasch
112! Initial revision
113!
114!
115! Description:
116! ------------
117!> Check control parameters and deduce further quantities.
118!------------------------------------------------------------------------------!
119 SUBROUTINE check_parameters
120
121
122    USE arrays_3d
123
124    USE basic_constants_and_equations_mod
125
126    USE bulk_cloud_model_mod,                                                  &
127        ONLY:  bulk_cloud_model
128
129    USE chem_modules
130
131    USE chemistry_model_mod,                                                   &
132        ONLY:  chem_boundary_conds
133
134    USE control_parameters
135
136    USE grid_variables
137
138    USE kinds
139
140    USE indices
141
142    USE model_1d_mod,                                                          &
143        ONLY:  damp_level_1d, damp_level_ind_1d
144
145    USE module_interface,                                                      &
146        ONLY:  module_interface_check_parameters,                              &
147               module_interface_check_data_output_ts,                          &
148               module_interface_check_data_output_pr,                          &
149               module_interface_check_data_output
150
151    USE netcdf_data_input_mod,                                                 &
152        ONLY:  init_model, input_pids_static, netcdf_data_input_check_dynamic, &
153               netcdf_data_input_check_static
154
155    USE netcdf_interface,                                                      &
156        ONLY:  dopr_unit, do2d_unit, do3d_unit, netcdf_data_format,            &
157               netcdf_data_format_string, dots_unit, heatflux_output_unit,     &
158               waterflux_output_unit, momentumflux_output_unit,                &
159               dots_max, dots_num, dots_label
160
161    USE particle_attributes,                                                   &
162        ONLY:  particle_advection, use_sgs_for_particles
163
164    USE pegrid
165
166    USE pmc_interface,                                                         &
167        ONLY:  cpl_id, nested_run
168
169    USE profil_parameter
170
171    USE statistics
172
173    USE subsidence_mod
174
175    USE transpose_indices
176
177#if defined( __parallel )
178    USE vertical_nesting_mod,                                                  &
179        ONLY:  vnested,                                                        &
180               vnest_check_parameters
181#endif
182
183
184    IMPLICIT NONE
185
186    CHARACTER (LEN=varnamelength)  ::  var           !< variable name
187    CHARACTER (LEN=7)   ::  unit                     !< unit of variable
188    CHARACTER (LEN=8)   ::  date                     !< current date string
189    CHARACTER (LEN=10)  ::  time                     !< current time string
190    CHARACTER (LEN=20)  ::  ensemble_string          !< string containing number of ensemble member
191    CHARACTER (LEN=15)  ::  nest_string              !< string containing id of nested domain
192    CHARACTER (LEN=40)  ::  coupling_string          !< string containing type of coupling
193    CHARACTER (LEN=100) ::  action                   !< flag string
194
195    INTEGER(iwp) ::  i                               !< loop index
196    INTEGER(iwp) ::  ilen                            !< string length
197    INTEGER(iwp) ::  j                               !< loop index
198    INTEGER(iwp) ::  k                               !< loop index
199    INTEGER(iwp) ::  kk                              !< loop index
200    INTEGER(iwp) ::  mid                             !< masked output running index
201    INTEGER(iwp) ::  netcdf_data_format_save         !< initial value of netcdf_data_format
202    INTEGER(iwp) ::  position                        !< index position of string
203
204    LOGICAL     ::  found                            !< flag, true if output variable is already marked for averaging
205
206    REAL(wp)    ::  gradient                         !< local gradient
207#if defined( __parallel )
208    REAL(wp)    ::  dt_spinup_max                    !< maximum spinup timestep in nested domains
209    REAL(wp)    ::  remote = 0.0_wp                  !< MPI id of remote processor
210    REAL(wp)    ::  spinup_time_max                  !< maximum spinup time in nested domains
211    REAL(wp)    ::  time_to_be_simulated_from_reference_point  !< time to be simulated from reference point
212#endif
213
214
215    CALL location_message( 'checking parameters', 'start' )
216!
217!-- At first, check static and dynamic input for consistency
218    CALL netcdf_data_input_check_dynamic
219    CALL netcdf_data_input_check_static
220!
221!-- Check for overlap combinations, which are not realized yet
222    IF ( transpose_compute_overlap  .AND. numprocs == 1 )  THEN
223          message_string = 'transpose-compute-overlap not implemented for single PE runs'
224          CALL message( 'check_parameters', 'PA0000', 1, 2, 0, 6, 0 )
225    ENDIF
226
227!
228!-- Check the coupling mode
229    IF ( coupling_mode /= 'uncoupled'            .AND.                         &
230         coupling_mode /= 'precursor_atmos'      .AND.                         &
231         coupling_mode /= 'precursor_ocean'      .AND.                         &
232         coupling_mode /= 'vnested_crse'         .AND.                         &
233         coupling_mode /= 'vnested_fine'         .AND.                         &
234         coupling_mode /= 'atmosphere_to_ocean'  .AND.                         &
235         coupling_mode /= 'ocean_to_atmosphere' )  THEN
236       message_string = 'illegal coupling mode: ' // TRIM( coupling_mode )
237       CALL message( 'check_parameters', 'PA0002', 1, 2, 0, 6, 0 )
238    ENDIF
239
240!
241!-- Check if humidity is set to TRUE in case of the atmospheric run (for coupled runs)
242    IF ( coupling_mode == 'atmosphere_to_ocean' .AND. .NOT. humidity) THEN
243       message_string = ' Humidity has to be set to .T. in the _p3d file ' //  &
244                        'for coupled runs between ocean and atmosphere.'
245       CALL message( 'check_parameters', 'PA0476', 1, 2, 0, 6, 0 )
246    ENDIF
247
248!
249!-- Check dt_coupling, restart_time, dt_restart, end_time, dx, dy, nx and ny
250    IF ( coupling_mode /= 'uncoupled'       .AND.                              &
251         coupling_mode(1:8) /= 'vnested_'   .AND.                              &
252         coupling_mode /= 'precursor_atmos' .AND.                              &
253         coupling_mode /= 'precursor_ocean' )  THEN
254
255       IF ( dt_coupling == 9999999.9_wp )  THEN
256          message_string = 'dt_coupling is not set but required for coup' //   &
257                           'ling mode "' //  TRIM( coupling_mode ) // '"'
258          CALL message( 'check_parameters', 'PA0003', 1, 2, 0, 6, 0 )
259       ENDIF
260
261#if defined( __parallel )
262
263
264       IF ( myid == 0 ) THEN
265          CALL MPI_SEND( dt_coupling, 1, MPI_REAL, target_id, 11, comm_inter,  &
266                         ierr )
267          CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 11, comm_inter,       &
268                         status, ierr )
269       ENDIF
270       CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr)
271
272       IF ( dt_coupling /= remote )  THEN
273          WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), &
274                 '": dt_coupling = ', dt_coupling, '& is not equal to ',       &
275                 'dt_coupling_remote = ', remote
276          CALL message( 'check_parameters', 'PA0004', 1, 2, 0, 6, 0 )
277       ENDIF
278       IF ( dt_coupling <= 0.0_wp )  THEN
279
280          IF ( myid == 0  ) THEN
281             CALL MPI_SEND( dt_max, 1, MPI_REAL, target_id, 19, comm_inter, ierr )
282             CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 19, comm_inter,    &
283                            status, ierr )
284          ENDIF
285          CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr)
286
287          dt_coupling = MAX( dt_max, remote )
288          WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), &
289                 '": dt_coupling <= 0.0 & is not allowed and is reset to ',    &
290                 'MAX(dt_max(A,O)) = ', dt_coupling
291          CALL message( 'check_parameters', 'PA0005', 0, 1, 0, 6, 0 )
292       ENDIF
293
294       IF ( myid == 0 ) THEN
295          CALL MPI_SEND( restart_time, 1, MPI_REAL, target_id, 12, comm_inter, &
296                         ierr )
297          CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 12, comm_inter,       &
298                         status, ierr )
299       ENDIF
300       CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr)
301
302       IF ( restart_time /= remote )  THEN
303          WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), &
304                 '": restart_time = ', restart_time, '& is not equal to ',     &
305                 'restart_time_remote = ', remote
306          CALL message( 'check_parameters', 'PA0006', 1, 2, 0, 6, 0 )
307       ENDIF
308
309       IF ( myid == 0 ) THEN
310          CALL MPI_SEND( dt_restart, 1, MPI_REAL, target_id, 13, comm_inter,   &
311                         ierr )
312          CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 13, comm_inter,       &
313                         status, ierr )
314       ENDIF
315       CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr)
316
317       IF ( dt_restart /= remote )  THEN
318          WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), &
319                 '": dt_restart = ', dt_restart, '& is not equal to ',         &
320                 'dt_restart_remote = ', remote
321          CALL message( 'check_parameters', 'PA0007', 1, 2, 0, 6, 0 )
322       ENDIF
323
324       time_to_be_simulated_from_reference_point = end_time-coupling_start_time
325
326       IF  ( myid == 0 ) THEN
327          CALL MPI_SEND( time_to_be_simulated_from_reference_point, 1,         &
328                         MPI_REAL, target_id, 14, comm_inter, ierr )
329          CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 14, comm_inter,       &
330                         status, ierr )
331       ENDIF
332       CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr)
333
334       IF ( time_to_be_simulated_from_reference_point /= remote )  THEN
335          WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), &
336                 '": time_to_be_simulated_from_reference_point = ',            &
337                 time_to_be_simulated_from_reference_point, '& is not equal ', &
338                 'to time_to_be_simulated_from_reference_point_remote = ',     &
339                 remote
340          CALL message( 'check_parameters', 'PA0008', 1, 2, 0, 6, 0 )
341       ENDIF
342
343       IF ( myid == 0 ) THEN
344          CALL MPI_SEND( dx, 1, MPI_REAL, target_id, 15, comm_inter, ierr )
345          CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 15, comm_inter,       &
346                                                             status, ierr )
347       ENDIF
348       CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr)
349
350
351       IF ( coupling_mode == 'atmosphere_to_ocean') THEN
352
353          IF ( dx < remote ) THEN
354             WRITE( message_string, * ) 'coupling mode "',                     &
355                   TRIM( coupling_mode ),                                      &
356           '": dx in Atmosphere is not equal to or not larger than dx in ocean'
357             CALL message( 'check_parameters', 'PA0009', 1, 2, 0, 6, 0 )
358          ENDIF
359
360          IF ( (nx_a+1)*dx /= (nx_o+1)*remote )  THEN
361             WRITE( message_string, * ) 'coupling mode "',                     &
362                    TRIM( coupling_mode ),                                     &
363             '": Domain size in x-direction is not equal in ocean and atmosphere'
364             CALL message( 'check_parameters', 'PA0010', 1, 2, 0, 6, 0 )
365          ENDIF
366
367       ENDIF
368
369       IF ( myid == 0) THEN
370          CALL MPI_SEND( dy, 1, MPI_REAL, target_id, 16, comm_inter, ierr )
371          CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 16, comm_inter,       &
372                         status, ierr )
373       ENDIF
374       CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr)
375
376       IF ( coupling_mode == 'atmosphere_to_ocean') THEN
377
378          IF ( dy < remote )  THEN
379             WRITE( message_string, * ) 'coupling mode "',                     &
380                    TRIM( coupling_mode ),                                     &
381                 '": dy in Atmosphere is not equal to or not larger than dy in ocean'
382             CALL message( 'check_parameters', 'PA0011', 1, 2, 0, 6, 0 )
383          ENDIF
384
385          IF ( (ny_a+1)*dy /= (ny_o+1)*remote )  THEN
386             WRITE( message_string, * ) 'coupling mode "',                     &
387                   TRIM( coupling_mode ),                                      &
388             '": Domain size in y-direction is not equal in ocean and atmosphere'
389             CALL message( 'check_parameters', 'PA0012', 1, 2, 0, 6, 0 )
390          ENDIF
391
392          IF ( MOD(nx_o+1,nx_a+1) /= 0 )  THEN
393             WRITE( message_string, * ) 'coupling mode "',                     &
394                   TRIM( coupling_mode ),                                      &
395             '": nx+1 in ocean is not divisible by nx+1 in',                   &
396             ' atmosphere without remainder'
397             CALL message( 'check_parameters', 'PA0339', 1, 2, 0, 6, 0 )
398          ENDIF
399
400          IF ( MOD(ny_o+1,ny_a+1) /= 0 )  THEN
401             WRITE( message_string, * ) 'coupling mode "',                     &
402                   TRIM( coupling_mode ),                                      &
403             '": ny+1 in ocean is not divisible by ny+1 in', &
404             ' atmosphere without remainder'
405
406             CALL message( 'check_parameters', 'PA0340', 1, 2, 0, 6, 0 )
407          ENDIF
408
409       ENDIF
410#else
411       WRITE( message_string, * ) 'coupling requires PALM to be compiled with',&
412            ' cpp-option "-D__parallel"'
413       CALL message( 'check_parameters', 'PA0141', 1, 2, 0, 6, 0 )
414#endif
415    ENDIF
416
417#if defined( __parallel )
418!
419!-- Exchange via intercommunicator
420    IF ( coupling_mode == 'atmosphere_to_ocean' .AND. myid == 0 )  THEN
421       CALL MPI_SEND( humidity, 1, MPI_LOGICAL, target_id, 19, comm_inter,     &
422                      ierr )
423    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' .AND. myid == 0)  THEN
424       CALL MPI_RECV( humidity_remote, 1, MPI_LOGICAL, target_id, 19,          &
425                      comm_inter, status, ierr )
426    ENDIF
427    CALL MPI_BCAST( humidity_remote, 1, MPI_LOGICAL, 0, comm2d, ierr)
428
429#endif
430
431!
432!-- User settings for restart times requires that "restart" has been given as
433!-- file activation string. Otherwise, binary output would not be saved by
434!-- palmrun.
435    IF (  ( restart_time /= 9999999.9_wp  .OR.  dt_restart /= 9999999.9_wp )   &
436         .AND.  .NOT. write_binary )  THEN
437       WRITE( message_string, * ) 'manual restart settings requires file ',    &
438                                  'activation string "restart"'
439       CALL message( 'check_parameters', 'PA0001', 1, 2, 0, 6, 0 )
440    ENDIF
441
442
443!
444!-- Generate the file header which is used as a header for most of PALM's
445!-- output files
446    CALL DATE_AND_TIME( date, time, run_zone )
447    run_date = date(1:4)//'-'//date(5:6)//'-'//date(7:8)
448    run_time = time(1:2)//':'//time(3:4)//':'//time(5:6)
449    IF ( coupling_mode == 'uncoupled' )  THEN
450       coupling_string = ''
451    ELSEIF ( coupling_mode == 'vnested_crse' )  THEN
452       coupling_string = ' nested (coarse)'
453    ELSEIF ( coupling_mode == 'vnested_fine' )  THEN
454       coupling_string = ' nested (fine)'
455    ELSEIF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
456       coupling_string = ' coupled (atmosphere)'
457    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
458       coupling_string = ' coupled (ocean)'
459    ENDIF
460    IF ( ensemble_member_nr /= 0 )  THEN
461       WRITE( ensemble_string, '(2X,A,I2.2)' )  'en-no: ', ensemble_member_nr
462    ELSE
463       ensemble_string = ''
464    ENDIF
465    IF ( nested_run )  THEN
466       WRITE( nest_string, '(2X,A,I2.2)' )  'nest-id: ', cpl_id
467    ELSE
468       nest_string = ''
469    ENDIF
470
471    WRITE ( run_description_header,                                            &
472            '(A,2X,A,2X,A,A,A,I2.2,A,A,A,2X,A,A,2X,A,1X,A)' )                  &
473          TRIM( version ), TRIM( revision ), 'run: ',                          &
474          TRIM( run_identifier ), '.', runnr, TRIM( coupling_string ),         &
475          TRIM( nest_string ), TRIM( ensemble_string), 'host: ', TRIM( host ), &
476          run_date, run_time
477
478!
479!-- Check the general loop optimization method
480    SELECT CASE ( TRIM( loop_optimization ) )
481
482       CASE ( 'cache', 'vector' )
483          CONTINUE
484
485       CASE DEFAULT
486          message_string = 'illegal value given for loop_optimization: "' //   &
487                           TRIM( loop_optimization ) // '"'
488          CALL message( 'check_parameters', 'PA0013', 1, 2, 0, 6, 0 )
489
490    END SELECT
491
492!
493!-- Check topography setting (check for illegal parameter combinations)
494    IF ( topography /= 'flat' )  THEN
495       action = ' '
496       IF ( scalar_advec /= 'pw-scheme' .AND. scalar_advec /= 'ws-scheme'      &
497          )  THEN
498          WRITE( action, '(A,A)' )  'scalar_advec = ', scalar_advec
499       ENDIF
500       IF ( momentum_advec /= 'pw-scheme' .AND. momentum_advec /= 'ws-scheme' )&
501       THEN
502          WRITE( action, '(A,A)' )  'momentum_advec = ', momentum_advec
503       ENDIF
504       IF ( psolver == 'sor' )  THEN
505          WRITE( action, '(A,A)' )  'psolver = ', psolver
506       ENDIF
507       IF ( sloping_surface )  THEN
508          WRITE( action, '(A)' )  'sloping surface = .TRUE.'
509       ENDIF
510       IF ( galilei_transformation )  THEN
511          WRITE( action, '(A)' )  'galilei_transformation = .TRUE.'
512       ENDIF
513       IF ( cloud_droplets )  THEN
514          WRITE( action, '(A)' )  'cloud_droplets = .TRUE.'
515       ENDIF
516       IF ( .NOT. constant_flux_layer .AND. topography /= 'closed_channel' )   &
517       THEN
518          WRITE( action, '(A)' )  'constant_flux_layer = .FALSE.'
519       ENDIF
520       IF ( action /= ' ' )  THEN
521          message_string = 'The specified topography does not allow ' //       &
522                           TRIM( action )
523          CALL message( 'check_parameters', 'PA0014', 1, 2, 0, 6, 0 )
524       ENDIF
525!
526!--    Check illegal/untested parameter combinations for closed channel
527       If ( topography == 'closed_channel' ) THEN
528          symmetry_flag = 1
529          message_string = 'Bottom and top boundary are treated equal'
530          CALL message( 'check_parameters', 'PA0699', 0, 0, 0, 6, 0 )
531
532          IF ( dz(1) /= dz(COUNT( dz /= -1.0_wp )) .OR.                        &
533               dz_stretch_level /= -9999999.9_wp) THEN
534             WRITE( message_string, * )  'dz should be equal close to the ' // &
535                                         'boundaries due to symmetrical problem'
536             CALL message( 'check_parameters', 'PA0700', 1, 2, 0, 6, 0 )
537          ENDIF
538
539          IF ( constant_flux_layer ) THEN
540             WRITE( message_string, * )  'A constant flux layer is not '//     &
541                                         'allowed if a closed channel '//      &
542                                         'shall be used'
543             CALL message( 'check_parameters', 'PA0701', 1, 2, 0, 6, 0 )
544          ENDIF
545
546          IF ( ocean_mode ) THEN
547             WRITE( message_string, * )  'The ocean mode is not allowed if '// &
548                                         'a closed channel shall be used'
549             CALL message( 'check_parameters', 'PA0702', 1, 2, 0, 6, 0 )
550          ENDIF
551
552          IF ( momentum_advec /= 'ws-scheme' .OR.                              &
553               scalar_advec /= 'ws-scheme' ) THEN
554             WRITE( message_string, * )  'A closed channel require the '//     &
555                                         'upwind scheme of Wicker and ' //     &
556                                         'Skamarock as the advection scheme'
557             CALL message( 'check_parameters', 'PA0703', 1, 2, 0, 6, 0 )
558          ENDIF
559       ENDIF
560    ENDIF
561
562!
563!-- Check approximation
564    IF ( TRIM( approximation ) /= 'boussinesq'   .AND.                         &
565         TRIM( approximation ) /= 'anelastic' )  THEN
566       message_string = 'unknown approximation: approximation = "' //          &
567                        TRIM( approximation ) // '"'
568       CALL message( 'check_parameters', 'PA0446', 1, 2, 0, 6, 0 )
569    ENDIF
570
571!
572!-- Check approximation requirements
573    IF ( TRIM( approximation ) == 'anelastic'   .AND.                          &
574         TRIM( momentum_advec ) /= 'ws-scheme' )  THEN
575       message_string = 'Anelastic approximation requires ' //                 &
576                        'momentum_advec = "ws-scheme"'
577       CALL message( 'check_parameters', 'PA0447', 1, 2, 0, 6, 0 )
578    ENDIF
579    IF ( TRIM( approximation ) == 'anelastic'   .AND.                          &
580         TRIM( psolver ) == 'multigrid' )  THEN
581       message_string = 'Anelastic approximation currently only supports ' //  &
582                        'psolver = "poisfft", ' //                             &
583                        'psolver = "sor" and ' //                              &
584                        'psolver = "multigrid_noopt"'
585       CALL message( 'check_parameters', 'PA0448', 1, 2, 0, 6, 0 )
586    ENDIF
587    IF ( TRIM( approximation ) == 'anelastic'   .AND.                          &
588         conserve_volume_flow )  THEN
589       message_string = 'Anelastic approximation is not allowed with ' //      &
590                        'conserve_volume_flow = .TRUE.'
591       CALL message( 'check_parameters', 'PA0449', 1, 2, 0, 6, 0 )
592    ENDIF
593
594!
595!-- Check flux input mode
596    IF ( TRIM( flux_input_mode ) /= 'dynamic'    .AND.                         &
597         TRIM( flux_input_mode ) /= 'kinematic'  .AND.                         &
598         TRIM( flux_input_mode ) /= 'approximation-specific' )  THEN
599       message_string = 'unknown flux input mode: flux_input_mode = "' //      &
600                        TRIM( flux_input_mode ) // '"'
601       CALL message( 'check_parameters', 'PA0450', 1, 2, 0, 6, 0 )
602    ENDIF
603!
604!-- Set flux input mode according to approximation if applicable
605    IF ( TRIM( flux_input_mode ) == 'approximation-specific' )  THEN
606       IF ( TRIM( approximation ) == 'anelastic' )  THEN
607          flux_input_mode = 'dynamic'
608       ELSEIF ( TRIM( approximation ) == 'boussinesq' )  THEN
609          flux_input_mode = 'kinematic'
610       ENDIF
611    ENDIF
612
613!
614!-- Check flux output mode
615    IF ( TRIM( flux_output_mode ) /= 'dynamic'    .AND.                        &
616         TRIM( flux_output_mode ) /= 'kinematic'  .AND.                        &
617         TRIM( flux_output_mode ) /= 'approximation-specific' )  THEN
618       message_string = 'unknown flux output mode: flux_output_mode = "' //    &
619                        TRIM( flux_output_mode ) // '"'
620       CALL message( 'check_parameters', 'PA0451', 1, 2, 0, 6, 0 )
621    ENDIF
622!
623!-- Set flux output mode according to approximation if applicable
624    IF ( TRIM( flux_output_mode ) == 'approximation-specific' )  THEN
625       IF ( TRIM( approximation ) == 'anelastic' )  THEN
626          flux_output_mode = 'dynamic'
627       ELSEIF ( TRIM( approximation ) == 'boussinesq' )  THEN
628          flux_output_mode = 'kinematic'
629       ENDIF
630    ENDIF
631
632
633!
634!-- When the land- or urban-surface model is used, the flux output must be
635!-- dynamic.
636    IF ( land_surface  .OR.  urban_surface )  THEN
637       flux_output_mode = 'dynamic'
638    ENDIF
639
640!
641!-- Set the flux output units according to flux_output_mode
642    IF ( TRIM( flux_output_mode ) == 'kinematic' ) THEN
643        heatflux_output_unit              = 'K m/s'
644        waterflux_output_unit             = 'kg/kg m/s'
645        momentumflux_output_unit          = 'm2/s2'
646    ELSEIF ( TRIM( flux_output_mode ) == 'dynamic' ) THEN
647        heatflux_output_unit              = 'W/m2'
648        waterflux_output_unit             = 'W/m2'
649        momentumflux_output_unit          = 'N/m2'
650    ENDIF
651
652!
653!-- set time series output units for fluxes
654    dots_unit(14:16) = TRIM( heatflux_output_unit )
655    dots_unit(21)    = TRIM( waterflux_output_unit )
656    dots_unit(19:20) = TRIM( momentumflux_output_unit )
657
658!
659!-- Add other module specific timeseries
660    CALL module_interface_check_data_output_ts( dots_max, dots_num, dots_label, dots_unit )
661
662!
663!-- Check if maximum number of allowed timeseries is exceeded
664    IF ( dots_num > dots_max )  THEN
665       WRITE( message_string, * ) 'number of time series quantities exceeds',  &
666                                  ' its maximum of dots_max = ', dots_max,     &
667                                  '&Please increase dots_max in modules.f90.'
668       CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 )
669    ENDIF
670
671!
672!-- Check whether there are any illegal values
673!-- Pressure solver:
674    IF ( psolver /= 'poisfft'  .AND.  psolver /= 'sor'  .AND.                  &
675         psolver /= 'multigrid'  .AND.  psolver /= 'multigrid_noopt' )  THEN
676       message_string = 'unknown solver for perturbation pressure: psolver' // &
677                        ' = "' // TRIM( psolver ) // '"'
678       CALL message( 'check_parameters', 'PA0016', 1, 2, 0, 6, 0 )
679    ENDIF
680
681    IF ( psolver(1:9) == 'multigrid' )  THEN
682       IF ( cycle_mg == 'w' )  THEN
683          gamma_mg = 2
684       ELSEIF ( cycle_mg == 'v' )  THEN
685          gamma_mg = 1
686       ELSE
687          message_string = 'unknown multigrid cycle: cycle_mg = "' //          &
688                           TRIM( cycle_mg ) // '"'
689          CALL message( 'check_parameters', 'PA0020', 1, 2, 0, 6, 0 )
690       ENDIF
691    ENDIF
692
693    IF ( fft_method /= 'singleton-algorithm'  .AND.                            &
694         fft_method /= 'temperton-algorithm'  .AND.                            &
695         fft_method /= 'fftw'                 .AND.                            &
696         fft_method /= 'system-specific' )  THEN
697       message_string = 'unknown fft-algorithm: fft_method = "' //             &
698                        TRIM( fft_method ) // '"'
699       CALL message( 'check_parameters', 'PA0021', 1, 2, 0, 6, 0 )
700    ENDIF
701
702    IF( momentum_advec == 'ws-scheme' .AND.                                    &
703        .NOT. call_psolver_at_all_substeps  ) THEN
704        message_string = 'psolver must be called at each RK3 substep when "'// &
705                      TRIM(momentum_advec) // ' "is used for momentum_advec'
706        CALL message( 'check_parameters', 'PA0344', 1, 2, 0, 6, 0 )
707    END IF
708!
709!-- Advection schemes:
710    IF ( momentum_advec /= 'pw-scheme'  .AND.                                  &
711         momentum_advec /= 'ws-scheme'  .AND.                                  &
712         momentum_advec /= 'up-scheme' )                                       &
713    THEN
714       message_string = 'unknown advection scheme: momentum_advec = "' //      &
715                        TRIM( momentum_advec ) // '"'
716       CALL message( 'check_parameters', 'PA0022', 1, 2, 0, 6, 0 )
717    ENDIF
718    IF ( ( momentum_advec == 'ws-scheme' .OR.  scalar_advec == 'ws-scheme' )   &
719           .AND. ( timestep_scheme == 'euler' .OR.                             &
720                   timestep_scheme == 'runge-kutta-2' ) )                      &
721    THEN
722       message_string = 'momentum_advec or scalar_advec = "'                   &
723         // TRIM( momentum_advec ) // '" is not allowed with ' //              &
724         'timestep_scheme = "' // TRIM( timestep_scheme ) // '"'
725       CALL message( 'check_parameters', 'PA0023', 1, 2, 0, 6, 0 )
726    ENDIF
727    IF ( scalar_advec /= 'pw-scheme'  .AND.  scalar_advec /= 'ws-scheme' .AND. &
728         scalar_advec /= 'bc-scheme' .AND. scalar_advec /= 'up-scheme' )       &
729    THEN
730       message_string = 'unknown advection scheme: scalar_advec = "' //        &
731                        TRIM( scalar_advec ) // '"'
732       CALL message( 'check_parameters', 'PA0024', 1, 2, 0, 6, 0 )
733    ENDIF
734    IF ( scalar_advec == 'bc-scheme'  .AND.  loop_optimization == 'cache' )    &
735    THEN
736       message_string = 'advection_scheme scalar_advec = "'                    &
737         // TRIM( scalar_advec ) // '" not implemented for ' //                &
738         'loop_optimization = "' // TRIM( loop_optimization ) // '"'
739       CALL message( 'check_parameters', 'PA0026', 1, 2, 0, 6, 0 )
740    ENDIF
741
742    IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets  .AND.             &
743         .NOT. use_upstream_for_tke  .AND.                                     &
744         scalar_advec /= 'ws-scheme'                                           &
745       )  THEN
746       use_upstream_for_tke = .TRUE.
747       message_string = 'use_upstream_for_tke is set to .TRUE. because ' //    &
748                        'use_sgs_for_particles = .TRUE. '          //          &
749                        'and scalar_advec /= ws-scheme'
750       CALL message( 'check_parameters', 'PA0025', 0, 1, 0, 6, 0 )
751    ENDIF
752
753!
754!-- Set LOGICAL switches to enhance performance
755    IF ( momentum_advec == 'ws-scheme' )  ws_scheme_mom = .TRUE.
756    IF ( scalar_advec   == 'ws-scheme' )  ws_scheme_sca = .TRUE.
757
758
759!
760!-- Timestep schemes:
761    SELECT CASE ( TRIM( timestep_scheme ) )
762
763       CASE ( 'euler' )
764          intermediate_timestep_count_max = 1
765
766       CASE ( 'runge-kutta-2' )
767          intermediate_timestep_count_max = 2
768
769       CASE ( 'runge-kutta-3' )
770          intermediate_timestep_count_max = 3
771
772       CASE DEFAULT
773          message_string = 'unknown timestep scheme: timestep_scheme = "' //   &
774                           TRIM( timestep_scheme ) // '"'
775          CALL message( 'check_parameters', 'PA0027', 1, 2, 0, 6, 0 )
776
777    END SELECT
778
779    IF ( (momentum_advec /= 'pw-scheme' .AND. momentum_advec /= 'ws-scheme')   &
780         .AND. timestep_scheme(1:5) == 'runge' ) THEN
781       message_string = 'momentum advection scheme "' // &
782                        TRIM( momentum_advec ) // '" & does not work with ' // &
783                        'timestep_scheme "' // TRIM( timestep_scheme ) // '"'
784       CALL message( 'check_parameters', 'PA0029', 1, 2, 0, 6, 0 )
785    ENDIF
786!
787!-- Check for proper settings for microphysics
788    IF ( bulk_cloud_model  .AND.  cloud_droplets )  THEN
789       message_string = 'bulk_cloud_model = .TRUE. is not allowed with ' //    &
790                        'cloud_droplets = .TRUE.'
791       CALL message( 'check_parameters', 'PA0442', 1, 2, 0, 6, 0 )
792    ENDIF
793
794!
795!-- Initializing actions must have been set by the user
796    IF ( TRIM( initializing_actions ) == '' )  THEN
797       message_string = 'no value specified for initializing_actions'
798       CALL message( 'check_parameters', 'PA0149', 1, 2, 0, 6, 0 )
799    ENDIF
800
801    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
802         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
803!
804!--    No restart run: several initialising actions are possible
805       action = initializing_actions
806       DO  WHILE ( TRIM( action ) /= '' )
807          position = INDEX( action, ' ' )
808          SELECT CASE ( action(1:position-1) )
809
810             CASE ( 'set_constant_profiles', 'set_1d-model_profiles',          &
811                    'by_user', 'initialize_vortex', 'initialize_ptanom',       &
812                    'initialize_bubble', 'inifor' )
813                action = action(position+1:)
814
815             CASE DEFAULT
816                message_string = 'initializing_action = "' //                  &
817                                 TRIM( action ) // '" unknown or not allowed'
818                CALL message( 'check_parameters', 'PA0030', 1, 2, 0, 6, 0 )
819
820          END SELECT
821       ENDDO
822    ENDIF
823
824    IF ( TRIM( initializing_actions ) == 'initialize_vortex'  .AND.            &
825         conserve_volume_flow ) THEN
826         message_string = 'initializing_actions = "initialize_vortex"' //      &
827                        ' is not allowed with conserve_volume_flow = .T.'
828       CALL message( 'check_parameters', 'PA0343', 1, 2, 0, 6, 0 )
829    ENDIF
830
831
832    IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0  .AND.    &
833         INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
834       message_string = 'initializing_actions = "set_constant_profiles"' //    &
835                        ' and "set_1d-model_profiles" are not allowed ' //     &
836                        'simultaneously'
837       CALL message( 'check_parameters', 'PA0031', 1, 2, 0, 6, 0 )
838    ENDIF
839
840    IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0  .AND.    &
841         INDEX( initializing_actions, 'by_user' ) /= 0 )  THEN
842       message_string = 'initializing_actions = "set_constant_profiles"' //    &
843                        ' and "by_user" are not allowed simultaneously'
844       CALL message( 'check_parameters', 'PA0032', 1, 2, 0, 6, 0 )
845    ENDIF
846
847    IF ( INDEX( initializing_actions, 'by_user' ) /= 0  .AND.                  &
848         INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
849       message_string = 'initializing_actions = "by_user" and ' //             &
850                        '"set_1d-model_profiles" are not allowed simultaneously'
851       CALL message( 'check_parameters', 'PA0033', 1, 2, 0, 6, 0 )
852    ENDIF
853!
854!-- In case of spinup and nested run, spinup end time must be identical
855!-- in order to have synchronously running simulations.
856    IF ( nested_run )  THEN
857#if defined( __parallel )
858       CALL MPI_ALLREDUCE( spinup_time, spinup_time_max, 1, MPI_REAL,          &
859                           MPI_MAX, MPI_COMM_WORLD, ierr )
860       CALL MPI_ALLREDUCE( dt_spinup,   dt_spinup_max,   1, MPI_REAL,          &
861                           MPI_MAX, MPI_COMM_WORLD, ierr )
862
863       IF ( spinup_time /= spinup_time_max  .OR.  dt_spinup /= dt_spinup_max ) &
864       THEN
865          message_string = 'In case of nesting, spinup_time and ' //           &
866                           'dt_spinup must be identical in all parent ' //     &
867                           'and child domains.'
868          CALL message( 'check_parameters', 'PA0489', 3, 2, 0, 6, 0 )
869       ENDIF
870#endif
871    ENDIF
872
873    IF ( bulk_cloud_model  .AND.  .NOT.  humidity )  THEN
874       WRITE( message_string, * ) 'bulk_cloud_model = ', bulk_cloud_model,     &
875              ' is not allowed with humidity = ', humidity
876       CALL message( 'check_parameters', 'PA0034', 1, 2, 0, 6, 0 )
877    ENDIF
878
879    IF ( humidity  .AND.  sloping_surface )  THEN
880       message_string = 'humidity = .TRUE. and sloping_surface = .TRUE. ' //   &
881                        'are not allowed simultaneously'
882       CALL message( 'check_parameters', 'PA0036', 1, 2, 0, 6, 0 )
883    ENDIF
884
885!-- Check the module settings
886    CALL module_interface_check_parameters
887
888!
889!-- In case of no restart run, check initialising parameters and calculate
890!-- further quantities
891    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
892
893!
894!--    Initial profiles for 1D and 3D model, respectively (u,v further below)
895       pt_init = pt_surface
896       IF ( humidity       )  q_init  = q_surface
897       IF ( passive_scalar )  s_init  = s_surface
898
899!--
900!--    If required, compute initial profile of the geostrophic wind
901!--    (component ug)
902       i = 1
903       gradient = 0.0_wp
904
905       IF ( .NOT. ocean_mode )  THEN
906
907          ug_vertical_gradient_level_ind(1) = 0
908          ug(0) = ug_surface
909          DO  k = 1, nzt+1
910             IF ( i < 11 )  THEN
911                IF ( ug_vertical_gradient_level(i) < zu(k)  .AND.              &
912                     ug_vertical_gradient_level(i) >= 0.0_wp )  THEN
913                   gradient = ug_vertical_gradient(i) / 100.0_wp
914                   ug_vertical_gradient_level_ind(i) = k - 1
915                   i = i + 1
916                ENDIF
917             ENDIF
918             IF ( gradient /= 0.0_wp )  THEN
919                IF ( k /= 1 )  THEN
920                   ug(k) = ug(k-1) + dzu(k) * gradient
921                ELSE
922                   ug(k) = ug_surface + dzu(k) * gradient
923                ENDIF
924             ELSE
925                ug(k) = ug(k-1)
926             ENDIF
927          ENDDO
928
929       ELSE
930
931          ug_vertical_gradient_level_ind(1) = nzt+1
932          ug(nzt+1) = ug_surface
933          DO  k = nzt, nzb, -1
934             IF ( i < 11 )  THEN
935                IF ( ug_vertical_gradient_level(i) > zu(k)  .AND.              &
936                     ug_vertical_gradient_level(i) <= 0.0_wp )  THEN
937                   gradient = ug_vertical_gradient(i) / 100.0_wp
938                   ug_vertical_gradient_level_ind(i) = k + 1
939                   i = i + 1
940                ENDIF
941             ENDIF
942             IF ( gradient /= 0.0_wp )  THEN
943                IF ( k /= nzt )  THEN
944                   ug(k) = ug(k+1) - dzu(k+1) * gradient
945                ELSE
946                   ug(k)   = ug_surface - 0.5_wp * dzu(k+1) * gradient
947                   ug(k+1) = ug_surface + 0.5_wp * dzu(k+1) * gradient
948                ENDIF
949             ELSE
950                ug(k) = ug(k+1)
951             ENDIF
952          ENDDO
953
954       ENDIF
955
956!
957!--    In case of no given gradients for ug, choose a zero gradient
958       IF ( ug_vertical_gradient_level(1) == -9999999.9_wp )  THEN
959          ug_vertical_gradient_level(1) = 0.0_wp
960       ENDIF
961
962!
963!--
964!--    If required, compute initial profile of the geostrophic wind
965!--    (component vg)
966       i = 1
967       gradient = 0.0_wp
968
969       IF ( .NOT. ocean_mode )  THEN
970
971          vg_vertical_gradient_level_ind(1) = 0
972          vg(0) = vg_surface
973          DO  k = 1, nzt+1
974             IF ( i < 11 )  THEN
975                IF ( vg_vertical_gradient_level(i) < zu(k)  .AND.              &
976                     vg_vertical_gradient_level(i) >= 0.0_wp )  THEN
977                   gradient = vg_vertical_gradient(i) / 100.0_wp
978                   vg_vertical_gradient_level_ind(i) = k - 1
979                   i = i + 1
980                ENDIF
981             ENDIF
982             IF ( gradient /= 0.0_wp )  THEN
983                IF ( k /= 1 )  THEN
984                   vg(k) = vg(k-1) + dzu(k) * gradient
985                ELSE
986                   vg(k) = vg_surface + dzu(k) * gradient
987                ENDIF
988             ELSE
989                vg(k) = vg(k-1)
990             ENDIF
991          ENDDO
992
993       ELSE
994
995          vg_vertical_gradient_level_ind(1) = nzt+1
996          vg(nzt+1) = vg_surface
997          DO  k = nzt, nzb, -1
998             IF ( i < 11 )  THEN
999                IF ( vg_vertical_gradient_level(i) > zu(k)  .AND.              &
1000                     vg_vertical_gradient_level(i) <= 0.0_wp )  THEN
1001                   gradient = vg_vertical_gradient(i) / 100.0_wp
1002                   vg_vertical_gradient_level_ind(i) = k + 1
1003                   i = i + 1
1004                ENDIF
1005             ENDIF
1006             IF ( gradient /= 0.0_wp )  THEN
1007                IF ( k /= nzt )  THEN
1008                   vg(k) = vg(k+1) - dzu(k+1) * gradient
1009                ELSE
1010                   vg(k)   = vg_surface - 0.5_wp * dzu(k+1) * gradient
1011                   vg(k+1) = vg_surface + 0.5_wp * dzu(k+1) * gradient
1012                ENDIF
1013             ELSE
1014                vg(k) = vg(k+1)
1015             ENDIF
1016          ENDDO
1017
1018       ENDIF
1019
1020!
1021!--    In case of no given gradients for vg, choose a zero gradient
1022       IF ( vg_vertical_gradient_level(1) == -9999999.9_wp )  THEN
1023          vg_vertical_gradient_level(1) = 0.0_wp
1024       ENDIF
1025
1026!
1027!--    Let the initial wind profiles be the calculated ug/vg profiles or
1028!--    interpolate them from wind profile data (if given)
1029       IF ( u_profile(1) == 9999999.9_wp  .AND.  v_profile(1) == 9999999.9_wp )  THEN
1030
1031          u_init = ug
1032          v_init = vg
1033
1034       ELSEIF ( u_profile(1) == 0.0_wp  .AND.  v_profile(1) == 0.0_wp )  THEN
1035
1036          IF ( uv_heights(1) /= 0.0_wp )  THEN
1037             message_string = 'uv_heights(1) must be 0.0'
1038             CALL message( 'check_parameters', 'PA0345', 1, 2, 0, 6, 0 )
1039          ENDIF
1040
1041          IF ( omega /= 0.0_wp )  THEN
1042             message_string = 'Coriolis force must be switched off (by setting omega=0.0)' //  &
1043                              ' when prescribing the forcing by u_profile and v_profile'
1044             CALL message( 'check_parameters', 'PA0347', 1, 2, 0, 6, 0 )
1045          ENDIF
1046
1047          use_prescribed_profile_data = .TRUE.
1048
1049          kk = 1
1050          u_init(0) = 0.0_wp
1051          v_init(0) = 0.0_wp
1052
1053          DO  k = 1, nz+1
1054
1055             IF ( kk < 200 )  THEN
1056                DO  WHILE ( uv_heights(kk+1) <= zu(k) )
1057                   kk = kk + 1
1058                   IF ( kk == 200 )  EXIT
1059                ENDDO
1060             ENDIF
1061
1062             IF ( kk < 200  .AND.  uv_heights(kk+1) /= 9999999.9_wp )  THEN
1063                u_init(k) = u_profile(kk) + ( zu(k) - uv_heights(kk) ) /       &
1064                                       ( uv_heights(kk+1) - uv_heights(kk) ) * &
1065                                       ( u_profile(kk+1) - u_profile(kk) )
1066                v_init(k) = v_profile(kk) + ( zu(k) - uv_heights(kk) ) /       &
1067                                       ( uv_heights(kk+1) - uv_heights(kk) ) * &
1068                                       ( v_profile(kk+1) - v_profile(kk) )
1069             ELSE
1070                u_init(k) = u_profile(kk)
1071                v_init(k) = v_profile(kk)
1072             ENDIF
1073
1074          ENDDO
1075
1076       ELSE
1077
1078          message_string = 'u_profile(1) and v_profile(1) must be 0.0'
1079          CALL message( 'check_parameters', 'PA0346', 1, 2, 0, 6, 0 )
1080
1081       ENDIF
1082
1083!
1084!--    Compute initial temperature profile using the given temperature gradients
1085       IF (  .NOT.  neutral )  THEN
1086          CALL init_vertical_profiles( pt_vertical_gradient_level_ind,          &
1087                                       pt_vertical_gradient_level,              &
1088                                       pt_vertical_gradient, pt_init,           &
1089                                       pt_surface, bc_pt_t_val )
1090       ENDIF
1091!
1092!--    Compute initial humidity profile using the given humidity gradients
1093       IF ( humidity )  THEN
1094          CALL init_vertical_profiles( q_vertical_gradient_level_ind,          &
1095                                       q_vertical_gradient_level,              &
1096                                       q_vertical_gradient, q_init,            &
1097                                       q_surface, bc_q_t_val )
1098       ENDIF
1099!
1100!--    Compute initial scalar profile using the given scalar gradients
1101       IF ( passive_scalar )  THEN
1102          CALL init_vertical_profiles( s_vertical_gradient_level_ind,          &
1103                                       s_vertical_gradient_level,              &
1104                                       s_vertical_gradient, s_init,            &
1105                                       s_surface, bc_s_t_val )
1106       ENDIF
1107!
1108!--    TODO
1109!--    Compute initial chemistry profile using the given chemical species gradients
1110!--    Russo: Is done in chem_init --> kanani: Revise
1111
1112    ENDIF
1113
1114!
1115!-- Check if the control parameter use_subsidence_tendencies is used correctly
1116    IF ( use_subsidence_tendencies  .AND.  .NOT.  large_scale_subsidence )  THEN
1117       message_string = 'The usage of use_subsidence_tendencies ' //           &
1118                            'requires large_scale_subsidence = .T..'
1119       CALL message( 'check_parameters', 'PA0396', 1, 2, 0, 6, 0 )
1120    ELSEIF ( use_subsidence_tendencies  .AND.  .NOT. large_scale_forcing )  THEN
1121       message_string = 'The usage of use_subsidence_tendencies ' //           &
1122                            'requires large_scale_forcing = .T..'
1123       CALL message( 'check_parameters', 'PA0397', 1, 2, 0, 6, 0 )
1124    ENDIF
1125
1126!
1127!-- Initialize large scale subsidence if required
1128    If ( large_scale_subsidence )  THEN
1129       IF ( subs_vertical_gradient_level(1) /= -9999999.9_wp  .AND.            &
1130                                     .NOT.  large_scale_forcing )  THEN
1131          CALL init_w_subsidence
1132       ENDIF
1133!
1134!--    In case large_scale_forcing is used, profiles for subsidence velocity
1135!--    are read in from file LSF_DATA
1136
1137       IF ( subs_vertical_gradient_level(1) == -9999999.9_wp  .AND.            &
1138            .NOT.  large_scale_forcing )  THEN
1139          message_string = 'There is no default large scale vertical ' //      &
1140                           'velocity profile set. Specify the subsidence ' //  &
1141                           'velocity profile via subs_vertical_gradient ' //   &
1142                           'and subs_vertical_gradient_level.'
1143          CALL message( 'check_parameters', 'PA0380', 1, 2, 0, 6, 0 )
1144       ENDIF
1145    ELSE
1146        IF ( subs_vertical_gradient_level(1) /= -9999999.9_wp )  THEN
1147           message_string = 'Enable usage of large scale subsidence by ' //    &
1148                            'setting large_scale_subsidence = .T..'
1149          CALL message( 'check_parameters', 'PA0381', 1, 2, 0, 6, 0 )
1150        ENDIF
1151    ENDIF
1152
1153!
1154!-- Overwrite parameters from namelist if necessary and compute Coriolis parameter.
1155!-- @todo - move initialization of f and fs to coriolis_mod.
1156    IF ( input_pids_static )  THEN
1157       latitude       = init_model%latitude
1158       longitude      = init_model%longitude
1159       rotation_angle = init_model%rotation_angle
1160    ENDIF
1161
1162    f  = 2.0_wp * omega * SIN( latitude / 180.0_wp * pi )
1163    fs = 2.0_wp * omega * COS( latitude / 180.0_wp * pi )
1164
1165!
1166!-- Check and set buoyancy related parameters and switches
1167    IF ( reference_state == 'horizontal_average' )  THEN
1168       CONTINUE
1169    ELSEIF ( reference_state == 'initial_profile' )  THEN
1170       use_initial_profile_as_reference = .TRUE.
1171    ELSEIF ( reference_state == 'single_value' )  THEN
1172       use_single_reference_value = .TRUE.
1173       IF ( pt_reference == 9999999.9_wp )  pt_reference = pt_surface
1174       vpt_reference = pt_reference * ( 1.0_wp + 0.61_wp * q_surface )
1175    ELSE
1176       message_string = 'illegal value for reference_state: "' //              &
1177                        TRIM( reference_state ) // '"'
1178       CALL message( 'check_parameters', 'PA0056', 1, 2, 0, 6, 0 )
1179    ENDIF
1180
1181!
1182!-- In case of a given slope, compute the relevant quantities
1183    IF ( alpha_surface /= 0.0_wp )  THEN
1184       IF ( ABS( alpha_surface ) > 90.0_wp )  THEN
1185          WRITE( message_string, * ) 'ABS( alpha_surface = ', alpha_surface,   &
1186                                     ' ) must be < 90.0'
1187          CALL message( 'check_parameters', 'PA0043', 1, 2, 0, 6, 0 )
1188       ENDIF
1189       sloping_surface = .TRUE.
1190       cos_alpha_surface = COS( alpha_surface / 180.0_wp * pi )
1191       sin_alpha_surface = SIN( alpha_surface / 180.0_wp * pi )
1192    ENDIF
1193
1194!
1195!-- Check time step and cfl_factor
1196    IF ( dt /= -1.0_wp )  THEN
1197       IF ( dt <= 0.0_wp )  THEN
1198          WRITE( message_string, * ) 'dt = ', dt , ' <= 0.0'
1199          CALL message( 'check_parameters', 'PA0044', 1, 2, 0, 6, 0 )
1200       ENDIF
1201       dt_3d = dt
1202       dt_fixed = .TRUE.
1203    ENDIF
1204
1205    IF ( cfl_factor <= 0.0_wp  .OR.  cfl_factor > 1.0_wp )  THEN
1206       IF ( cfl_factor == -1.0_wp )  THEN
1207          IF ( timestep_scheme == 'runge-kutta-2' )  THEN
1208             cfl_factor = 0.8_wp
1209          ELSEIF ( timestep_scheme == 'runge-kutta-3' )  THEN
1210             cfl_factor = 0.9_wp
1211          ELSE
1212             cfl_factor = 0.9_wp
1213          ENDIF
1214       ELSE
1215          WRITE( message_string, * ) 'cfl_factor = ', cfl_factor,              &
1216                 ' out of range &0.0 < cfl_factor <= 1.0 is required'
1217          CALL message( 'check_parameters', 'PA0045', 1, 2, 0, 6, 0 )
1218       ENDIF
1219    ENDIF
1220
1221!
1222!-- Store simulated time at begin
1223    simulated_time_at_begin = simulated_time
1224
1225!
1226!-- Store reference time for coupled runs and change the coupling flag,
1227!-- if ...
1228    IF ( simulated_time == 0.0_wp )  THEN
1229       IF ( coupling_start_time == 0.0_wp )  THEN
1230          time_since_reference_point = 0.0_wp
1231       ELSEIF ( time_since_reference_point < 0.0_wp )  THEN
1232          run_coupled = .FALSE.
1233       ENDIF
1234    ENDIF
1235
1236!
1237!-- Set wind speed in the Galilei-transformed system
1238    IF ( galilei_transformation )  THEN
1239       IF ( use_ug_for_galilei_tr                    .AND.                     &
1240            ug_vertical_gradient_level(1) == 0.0_wp  .AND.                     &
1241            ug_vertical_gradient(1) == 0.0_wp        .AND.                     &
1242            vg_vertical_gradient_level(1) == 0.0_wp  .AND.                     &
1243            vg_vertical_gradient(1) == 0.0_wp )  THEN
1244          u_gtrans = ug_surface * 0.6_wp
1245          v_gtrans = vg_surface * 0.6_wp
1246       ELSEIF ( use_ug_for_galilei_tr  .AND.                                   &
1247                ( ug_vertical_gradient_level(1) /= 0.0_wp  .OR.                &
1248                ug_vertical_gradient(1) /= 0.0_wp ) )  THEN
1249          message_string = 'baroclinity (ug) not allowed simultaneously' //    &
1250                           ' with galilei transformation'
1251          CALL message( 'check_parameters', 'PA0046', 1, 2, 0, 6, 0 )
1252       ELSEIF ( use_ug_for_galilei_tr  .AND.                                   &
1253                ( vg_vertical_gradient_level(1) /= 0.0_wp  .OR.                &
1254                vg_vertical_gradient(1) /= 0.0_wp ) )  THEN
1255          message_string = 'baroclinity (vg) not allowed simultaneously' //    &
1256                           ' with galilei transformation'
1257          CALL message( 'check_parameters', 'PA0047', 1, 2, 0, 6, 0 )
1258       ELSE
1259          message_string = 'variable translation speed used for Galilei-' //   &
1260             'transformation, which may cause & instabilities in stably ' //   &
1261             'stratified regions'
1262          CALL message( 'check_parameters', 'PA0048', 0, 1, 0, 6, 0 )
1263       ENDIF
1264    ENDIF
1265
1266!
1267!-- In case of using a prandtl-layer, calculated (or prescribed) surface
1268!-- fluxes have to be used in the diffusion-terms
1269    IF ( constant_flux_layer )  use_surface_fluxes = .TRUE.
1270
1271!
1272!-- Check boundary conditions and set internal variables:
1273!-- Attention: the lateral boundary conditions have been already checked in
1274!-- parin
1275!
1276!-- Non-cyclic lateral boundaries require the multigrid method and Piascek-
1277!-- Willimas or Wicker - Skamarock advection scheme. Several schemes
1278!-- and tools do not work with non-cyclic boundary conditions.
1279    IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
1280       IF ( psolver(1:9) /= 'multigrid' )  THEN
1281          message_string = 'non-cyclic lateral boundaries do not allow ' //    &
1282                           'psolver = "' // TRIM( psolver ) // '"'
1283          CALL message( 'check_parameters', 'PA0051', 1, 2, 0, 6, 0 )
1284       ENDIF
1285       IF ( momentum_advec /= 'pw-scheme'  .AND.                               &
1286            momentum_advec /= 'ws-scheme' )  THEN
1287
1288          message_string = 'non-cyclic lateral boundaries do not allow ' //    &
1289                           'momentum_advec = "' // TRIM( momentum_advec ) // '"'
1290          CALL message( 'check_parameters', 'PA0052', 1, 2, 0, 6, 0 )
1291       ENDIF
1292       IF ( scalar_advec /= 'pw-scheme'  .AND.                                 &
1293            scalar_advec /= 'ws-scheme' )  THEN
1294          message_string = 'non-cyclic lateral boundaries do not allow ' //    &
1295                           'scalar_advec = "' // TRIM( scalar_advec ) // '"'
1296          CALL message( 'check_parameters', 'PA0053', 1, 2, 0, 6, 0 )
1297       ENDIF
1298       IF ( galilei_transformation )  THEN
1299          message_string = 'non-cyclic lateral boundaries do not allow ' //    &
1300                           'galilei_transformation = .T.'
1301          CALL message( 'check_parameters', 'PA0054', 1, 2, 0, 6, 0 )
1302       ENDIF
1303    ENDIF
1304
1305!
1306!-- Bottom boundary condition for the turbulent Kinetic energy
1307    IF ( bc_e_b == 'neumann' )  THEN
1308       ibc_e_b = 1
1309    ELSEIF ( bc_e_b == '(u*)**2+neumann' )  THEN
1310       ibc_e_b = 2
1311       IF ( .NOT. constant_flux_layer )  THEN
1312          bc_e_b = 'neumann'
1313          ibc_e_b = 1
1314          message_string = 'boundary condition bc_e_b changed to "' //         &
1315                           TRIM( bc_e_b ) // '"'
1316          CALL message( 'check_parameters', 'PA0057', 0, 1, 0, 6, 0 )
1317       ENDIF
1318    ELSE
1319       message_string = 'unknown boundary condition: bc_e_b = "' //            &
1320                        TRIM( bc_e_b ) // '"'
1321       CALL message( 'check_parameters', 'PA0058', 1, 2, 0, 6, 0 )
1322    ENDIF
1323
1324!
1325!-- Boundary conditions for perturbation pressure
1326    IF ( bc_p_b == 'dirichlet' )  THEN
1327       ibc_p_b = 0
1328    ELSEIF ( bc_p_b == 'neumann' )  THEN
1329       ibc_p_b = 1
1330    ELSE
1331       message_string = 'unknown boundary condition: bc_p_b = "' //            &
1332                        TRIM( bc_p_b ) // '"'
1333       CALL message( 'check_parameters', 'PA0059', 1, 2, 0, 6, 0 )
1334    ENDIF
1335
1336    IF ( bc_p_t == 'dirichlet' )  THEN
1337       ibc_p_t = 0
1338!-- TO_DO: later set bc_p_t to neumann before, in case of nested domain
1339    ELSEIF ( bc_p_t == 'neumann' .OR. bc_p_t == 'nested' )  THEN
1340       ibc_p_t = 1
1341    ELSE
1342       message_string = 'unknown boundary condition: bc_p_t = "' //            &
1343                        TRIM( bc_p_t ) // '"'
1344       CALL message( 'check_parameters', 'PA0061', 1, 2, 0, 6, 0 )
1345    ENDIF
1346
1347!
1348!-- Boundary conditions for potential temperature
1349    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
1350       ibc_pt_b = 2
1351    ELSE
1352       IF ( bc_pt_b == 'dirichlet' )  THEN
1353          ibc_pt_b = 0
1354       ELSEIF ( bc_pt_b == 'neumann' )  THEN
1355          ibc_pt_b = 1
1356       ELSE
1357          message_string = 'unknown boundary condition: bc_pt_b = "' //        &
1358                           TRIM( bc_pt_b ) // '"'
1359          CALL message( 'check_parameters', 'PA0062', 1, 2, 0, 6, 0 )
1360       ENDIF
1361    ENDIF
1362
1363    IF ( bc_pt_t == 'dirichlet' )  THEN
1364       ibc_pt_t = 0
1365    ELSEIF ( bc_pt_t == 'neumann' )  THEN
1366       ibc_pt_t = 1
1367    ELSEIF ( bc_pt_t == 'initial_gradient' )  THEN
1368       ibc_pt_t = 2
1369    ELSEIF ( bc_pt_t == 'nested'  .OR.  bc_pt_t == 'nesting_offline' )  THEN
1370       ibc_pt_t = 3
1371    ELSE
1372       message_string = 'unknown boundary condition: bc_pt_t = "' //           &
1373                        TRIM( bc_pt_t ) // '"'
1374       CALL message( 'check_parameters', 'PA0063', 1, 2, 0, 6, 0 )
1375    ENDIF
1376
1377    IF ( ANY( wall_heatflux /= 0.0_wp )  .AND.                        &
1378         surface_heatflux == 9999999.9_wp )  THEN
1379       message_string = 'wall_heatflux additionally requires ' //     &
1380                        'setting of surface_heatflux'
1381       CALL message( 'check_parameters', 'PA0443', 1, 2, 0, 6, 0 )
1382    ENDIF
1383
1384!
1385!   This IF clause needs revision, got too complex!!
1386    IF ( surface_heatflux == 9999999.9_wp  )  THEN
1387       constant_heatflux = .FALSE.
1388       IF ( large_scale_forcing  .OR.  land_surface  .OR.  urban_surface )  THEN
1389          IF ( ibc_pt_b == 0 )  THEN
1390             constant_heatflux = .FALSE.
1391          ELSEIF ( ibc_pt_b == 1 )  THEN
1392             constant_heatflux = .TRUE.
1393             surface_heatflux = 0.0_wp
1394          ENDIF
1395       ENDIF
1396    ELSE
1397       constant_heatflux = .TRUE.
1398    ENDIF
1399
1400    IF ( top_heatflux     == 9999999.9_wp )  constant_top_heatflux = .FALSE.
1401
1402    IF ( neutral )  THEN
1403
1404       IF ( surface_heatflux /= 0.0_wp  .AND.                                  &
1405            surface_heatflux /= 9999999.9_wp )  THEN
1406          message_string = 'heatflux must not be set for pure neutral flow'
1407          CALL message( 'check_parameters', 'PA0351', 1, 2, 0, 6, 0 )
1408       ENDIF
1409
1410       IF ( top_heatflux /= 0.0_wp  .AND.  top_heatflux /= 9999999.9_wp )      &
1411       THEN
1412          message_string = 'heatflux must not be set for pure neutral flow'
1413          CALL message( 'check_parameters', 'PA0351', 1, 2, 0, 6, 0 )
1414       ENDIF
1415
1416    ENDIF
1417
1418    IF ( top_momentumflux_u /= 9999999.9_wp  .AND.                             &
1419         top_momentumflux_v /= 9999999.9_wp )  THEN
1420       constant_top_momentumflux = .TRUE.
1421    ELSEIF (  .NOT. ( top_momentumflux_u == 9999999.9_wp  .AND.                &
1422           top_momentumflux_v == 9999999.9_wp ) )  THEN
1423       message_string = 'both, top_momentumflux_u AND top_momentumflux_v ' //  &
1424                        'must be set'
1425       CALL message( 'check_parameters', 'PA0064', 1, 2, 0, 6, 0 )
1426    ENDIF
1427
1428!
1429!-- A given surface temperature implies Dirichlet boundary condition for
1430!-- temperature. In this case specification of a constant heat flux is
1431!-- forbidden.
1432    IF ( ibc_pt_b == 0  .AND.  constant_heatflux  .AND.                        &
1433         surface_heatflux /= 0.0_wp )  THEN
1434       message_string = 'boundary_condition: bc_pt_b = "' // TRIM( bc_pt_b ) //&
1435                        '& is not allowed with constant_heatflux = .TRUE.'
1436       CALL message( 'check_parameters', 'PA0065', 1, 2, 0, 6, 0 )
1437    ENDIF
1438    IF ( constant_heatflux  .AND.  pt_surface_initial_change /= 0.0_wp )  THEN
1439       WRITE ( message_string, * )  'constant_heatflux = .TRUE. is not allo',  &
1440               'wed with pt_surface_initial_change (/=0) = ',                  &
1441               pt_surface_initial_change
1442       CALL message( 'check_parameters', 'PA0066', 1, 2, 0, 6, 0 )
1443    ENDIF
1444
1445!
1446!-- A given temperature at the top implies Dirichlet boundary condition for
1447!-- temperature. In this case specification of a constant heat flux is
1448!-- forbidden.
1449    IF ( ibc_pt_t == 0  .AND.  constant_top_heatflux  .AND.                    &
1450         top_heatflux /= 0.0_wp )  THEN
1451       message_string = 'boundary_condition: bc_pt_t = "' // TRIM( bc_pt_t ) //&
1452                        '" is not allowed with constant_top_heatflux = .TRUE.'
1453       CALL message( 'check_parameters', 'PA0067', 1, 2, 0, 6, 0 )
1454    ENDIF
1455
1456!
1457!-- Set boundary conditions for total water content
1458    IF ( humidity )  THEN
1459
1460       IF ( ANY( wall_humidityflux /= 0.0_wp )  .AND.                        &
1461            surface_waterflux == 9999999.9_wp )  THEN
1462          message_string = 'wall_humidityflux additionally requires ' //     &
1463                           'setting of surface_waterflux'
1464          CALL message( 'check_parameters', 'PA0444', 1, 2, 0, 6, 0 )
1465       ENDIF
1466
1467       CALL set_bc_scalars( 'q', bc_q_b, bc_q_t, ibc_q_b, ibc_q_t,           &
1468                            'PA0071', 'PA0072' )
1469
1470       IF ( surface_waterflux == 9999999.9_wp  )  THEN
1471          constant_waterflux = .FALSE.
1472          IF ( large_scale_forcing .OR. land_surface )  THEN
1473             IF ( ibc_q_b == 0 )  THEN
1474                constant_waterflux = .FALSE.
1475             ELSEIF ( ibc_q_b == 1 )  THEN
1476                constant_waterflux = .TRUE.
1477             ENDIF
1478          ENDIF
1479       ELSE
1480          constant_waterflux = .TRUE.
1481       ENDIF
1482
1483       CALL check_bc_scalars( 'q', bc_q_b, ibc_q_b, 'PA0073', 'PA0074',        &
1484                              constant_waterflux, q_surface_initial_change )
1485
1486    ENDIF
1487
1488    IF ( passive_scalar )  THEN
1489
1490       IF ( ANY( wall_scalarflux /= 0.0_wp )  .AND.                            &
1491            surface_scalarflux == 9999999.9_wp )  THEN
1492          message_string = 'wall_scalarflux additionally requires ' //         &
1493                           'setting of surface_scalarflux'
1494          CALL message( 'check_parameters', 'PA0445', 1, 2, 0, 6, 0 )
1495       ENDIF
1496
1497       IF ( surface_scalarflux == 9999999.9_wp )  constant_scalarflux = .FALSE.
1498
1499       CALL set_bc_scalars( 's', bc_s_b, bc_s_t, ibc_s_b, ibc_s_t,             &
1500                            'PA0071', 'PA0072' )
1501
1502       CALL check_bc_scalars( 's', bc_s_b, ibc_s_b, 'PA0073', 'PA0074',        &
1503                              constant_scalarflux, s_surface_initial_change )
1504
1505       IF ( top_scalarflux == 9999999.9_wp )  constant_top_scalarflux = .FALSE.
1506!
1507!--    A fixed scalar concentration at the top implies Dirichlet boundary
1508!--    condition for scalar. Hence, in this case specification of a constant
1509!--    scalar flux is forbidden.
1510       IF ( ( ibc_s_t == 0 .OR. ibc_s_t == 2 )  .AND.  constant_top_scalarflux &
1511               .AND.  top_scalarflux /= 0.0_wp )  THEN
1512          message_string = 'boundary condition: bc_s_t = "' //                 &
1513                           TRIM( bc_s_t ) // '" is not allowed with ' //       &
1514                           'top_scalarflux /= 0.0'
1515          CALL message( 'check_parameters', 'PA0441', 1, 2, 0, 6, 0 )
1516       ENDIF
1517    ENDIF
1518
1519!
1520!-- Boundary conditions for chemical species
1521    IF ( air_chemistry )  CALL chem_boundary_conds( 'init' )
1522
1523!
1524!-- Boundary conditions for horizontal components of wind speed
1525    IF ( bc_uv_b == 'dirichlet' )  THEN
1526       ibc_uv_b = 0
1527    ELSEIF ( bc_uv_b == 'neumann' )  THEN
1528       ibc_uv_b = 1
1529       IF ( constant_flux_layer )  THEN
1530          message_string = 'boundary condition: bc_uv_b = "' //                &
1531               TRIM( bc_uv_b ) // '" is not allowed with constant_flux_layer'  &
1532               // ' = .TRUE.'
1533          CALL message( 'check_parameters', 'PA0075', 1, 2, 0, 6, 0 )
1534       ENDIF
1535    ELSE
1536       message_string = 'unknown boundary condition: bc_uv_b = "' //           &
1537                        TRIM( bc_uv_b ) // '"'
1538       CALL message( 'check_parameters', 'PA0076', 1, 2, 0, 6, 0 )
1539    ENDIF
1540!
1541!-- In case of coupled simulations u and v at the ground in atmosphere will be
1542!-- assigned with the u and v values of the ocean surface
1543    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
1544       ibc_uv_b = 2
1545    ENDIF
1546
1547    IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1548       bc_uv_t = 'neumann'
1549       ibc_uv_t = 1
1550    ELSE
1551       IF ( bc_uv_t == 'dirichlet' .OR. bc_uv_t == 'dirichlet_0' )  THEN
1552          ibc_uv_t = 0
1553          IF ( bc_uv_t == 'dirichlet_0' )  THEN
1554!
1555!--          Velocities for the initial u,v-profiles are set zero at the top
1556!--          in case of dirichlet_0 conditions
1557             u_init(nzt+1)    = 0.0_wp
1558             v_init(nzt+1)    = 0.0_wp
1559          ENDIF
1560       ELSEIF ( bc_uv_t == 'neumann' )  THEN
1561          ibc_uv_t = 1
1562       ELSEIF ( bc_uv_t == 'nested'  .OR.  bc_uv_t == 'nesting_offline' )  THEN
1563          ibc_uv_t = 3
1564       ELSE
1565          message_string = 'unknown boundary condition: bc_uv_t = "' //        &
1566                           TRIM( bc_uv_t ) // '"'
1567          CALL message( 'check_parameters', 'PA0077', 1, 2, 0, 6, 0 )
1568       ENDIF
1569    ENDIF
1570
1571!
1572!-- Compute and check, respectively, the Rayleigh Damping parameter
1573    IF ( rayleigh_damping_factor == -1.0_wp )  THEN
1574       rayleigh_damping_factor = 0.0_wp
1575    ELSE
1576       IF ( rayleigh_damping_factor < 0.0_wp  .OR.                             &
1577            rayleigh_damping_factor > 1.0_wp )  THEN
1578          WRITE( message_string, * )  'rayleigh_damping_factor = ',            &
1579                              rayleigh_damping_factor, ' out of range [0.0,1.0]'
1580          CALL message( 'check_parameters', 'PA0078', 1, 2, 0, 6, 0 )
1581       ENDIF
1582    ENDIF
1583
1584    IF ( rayleigh_damping_height == -1.0_wp )  THEN
1585       IF (  .NOT.  ocean_mode )  THEN
1586          rayleigh_damping_height = 0.66666666666_wp * zu(nzt)
1587       ELSE
1588          rayleigh_damping_height = 0.66666666666_wp * zu(nzb)
1589       ENDIF
1590    ELSE
1591       IF (  .NOT.  ocean_mode )  THEN
1592          IF ( rayleigh_damping_height < 0.0_wp  .OR.                          &
1593               rayleigh_damping_height > zu(nzt) )  THEN
1594             WRITE( message_string, * )  'rayleigh_damping_height = ',         &
1595                   rayleigh_damping_height, ' out of range [0.0,', zu(nzt), ']'
1596             CALL message( 'check_parameters', 'PA0079', 1, 2, 0, 6, 0 )
1597          ENDIF
1598       ELSE
1599          IF ( rayleigh_damping_height > 0.0_wp  .OR.                          &
1600               rayleigh_damping_height < zu(nzb) )  THEN
1601             WRITE( message_string, * )  'rayleigh_damping_height = ',         &
1602                   rayleigh_damping_height, ' out of range [0.0,', zu(nzb), ']'
1603             CALL message( 'check_parameters', 'PA0079', 1, 2, 0, 6, 0 )
1604          ENDIF
1605       ENDIF
1606    ENDIF
1607
1608!
1609!-- Check number of chosen statistic regions
1610    IF ( statistic_regions < 0 )  THEN
1611       WRITE ( message_string, * ) 'number of statistic_regions = ',           &
1612                   statistic_regions+1, ' is not allowed'
1613       CALL message( 'check_parameters', 'PA0082', 1, 2, 0, 6, 0 )
1614    ENDIF
1615    IF ( normalizing_region > statistic_regions  .OR.                          &
1616         normalizing_region < 0)  THEN
1617       WRITE ( message_string, * ) 'normalizing_region = ',                    &
1618                normalizing_region, ' must be >= 0 and <= ',statistic_regions, &
1619                ' (value of statistic_regions)'
1620       CALL message( 'check_parameters', 'PA0083', 1, 2, 0, 6, 0 )
1621    ENDIF
1622
1623!
1624!-- Set the default intervals for data output, if necessary
1625!-- NOTE: dt_dosp has already been set in spectra_parin
1626    IF ( dt_data_output /= 9999999.9_wp )  THEN
1627       IF ( dt_dopr           == 9999999.9_wp )  dt_dopr           = dt_data_output
1628       IF ( dt_dopts          == 9999999.9_wp )  dt_dopts          = dt_data_output
1629       IF ( dt_do2d_xy        == 9999999.9_wp )  dt_do2d_xy        = dt_data_output
1630       IF ( dt_do2d_xz        == 9999999.9_wp )  dt_do2d_xz        = dt_data_output
1631       IF ( dt_do2d_yz        == 9999999.9_wp )  dt_do2d_yz        = dt_data_output
1632       IF ( dt_do3d           == 9999999.9_wp )  dt_do3d           = dt_data_output
1633       IF ( dt_data_output_av == 9999999.9_wp )  dt_data_output_av = dt_data_output
1634       DO  mid = 1, max_masks
1635          IF ( dt_domask(mid) == 9999999.9_wp )  dt_domask(mid)    = dt_data_output
1636       ENDDO
1637    ENDIF
1638
1639!
1640!-- Set the default skip time intervals for data output, if necessary
1641    IF ( skip_time_dopr    == 9999999.9_wp )                                   &
1642                                       skip_time_dopr    = skip_time_data_output
1643    IF ( skip_time_do2d_xy == 9999999.9_wp )                                   &
1644                                       skip_time_do2d_xy = skip_time_data_output
1645    IF ( skip_time_do2d_xz == 9999999.9_wp )                                   &
1646                                       skip_time_do2d_xz = skip_time_data_output
1647    IF ( skip_time_do2d_yz == 9999999.9_wp )                                   &
1648                                       skip_time_do2d_yz = skip_time_data_output
1649    IF ( skip_time_do3d    == 9999999.9_wp )                                   &
1650                                       skip_time_do3d    = skip_time_data_output
1651    IF ( skip_time_data_output_av == 9999999.9_wp )                            &
1652                                skip_time_data_output_av = skip_time_data_output
1653    DO  mid = 1, max_masks
1654       IF ( skip_time_domask(mid) == 9999999.9_wp )                            &
1655                                skip_time_domask(mid)    = skip_time_data_output
1656    ENDDO
1657
1658!
1659!-- Check the average intervals (first for 3d-data, then for profiles)
1660    IF ( averaging_interval > dt_data_output_av )  THEN
1661       WRITE( message_string, * )  'averaging_interval = ',                    &
1662             averaging_interval, ' must be <= dt_data_output_av = ',           &
1663             dt_data_output_av
1664       CALL message( 'check_parameters', 'PA0085', 1, 2, 0, 6, 0 )
1665    ENDIF
1666
1667    IF ( averaging_interval_pr == 9999999.9_wp )  THEN
1668       averaging_interval_pr = averaging_interval
1669    ENDIF
1670
1671    IF ( averaging_interval_pr > dt_dopr )  THEN
1672       WRITE( message_string, * )  'averaging_interval_pr = ',                 &
1673             averaging_interval_pr, ' must be <= dt_dopr = ', dt_dopr
1674       CALL message( 'check_parameters', 'PA0086', 1, 2, 0, 6, 0 )
1675    ENDIF
1676
1677!
1678!-- Set the default interval for profiles entering the temporal average
1679    IF ( dt_averaging_input_pr == 9999999.9_wp )  THEN
1680       dt_averaging_input_pr = dt_averaging_input
1681    ENDIF
1682
1683!
1684!-- Set the default interval for the output of timeseries to a reasonable
1685!-- value (tries to minimize the number of calls of flow_statistics)
1686    IF ( dt_dots == 9999999.9_wp )  THEN
1687       IF ( averaging_interval_pr == 0.0_wp )  THEN
1688          dt_dots = MIN( dt_run_control, dt_dopr )
1689       ELSE
1690          dt_dots = MIN( dt_run_control, dt_averaging_input_pr )
1691       ENDIF
1692    ENDIF
1693
1694!
1695!-- Check the sample rate for averaging (first for 3d-data, then for profiles)
1696    IF ( dt_averaging_input > averaging_interval )  THEN
1697       WRITE( message_string, * )  'dt_averaging_input = ',                    &
1698                dt_averaging_input, ' must be <= averaging_interval = ',       &
1699                averaging_interval
1700       CALL message( 'check_parameters', 'PA0088', 1, 2, 0, 6, 0 )
1701    ENDIF
1702
1703    IF ( dt_averaging_input_pr > averaging_interval_pr )  THEN
1704       WRITE( message_string, * )  'dt_averaging_input_pr = ',                 &
1705                dt_averaging_input_pr, ' must be <= averaging_interval_pr = ', &
1706                averaging_interval_pr
1707       CALL message( 'check_parameters', 'PA0089', 1, 2, 0, 6, 0 )
1708    ENDIF
1709
1710!
1711!-- Determine the number of output profiles and check whether they are
1712!-- permissible
1713    DO  WHILE ( data_output_pr(dopr_n+1) /= '          ' )
1714
1715       dopr_n = dopr_n + 1
1716       i = dopr_n
1717
1718!
1719!--    Determine internal profile number (for hom, homs)
1720!--    and store height levels
1721       SELECT CASE ( TRIM( data_output_pr(i) ) )
1722
1723          CASE ( 'u', '#u' )
1724             dopr_index(i) = 1
1725             dopr_unit(i)  = 'm/s'
1726             hom(:,2,1,:)  = SPREAD( zu, 2, statistic_regions+1 )
1727             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1728                dopr_initial_index(i) = 5
1729                hom(:,2,5,:)          = SPREAD( zu, 2, statistic_regions+1 )
1730                data_output_pr(i)     = data_output_pr(i)(2:)
1731             ENDIF
1732
1733          CASE ( 'v', '#v' )
1734             dopr_index(i) = 2
1735             dopr_unit(i)  = 'm/s'
1736             hom(:,2,2,:)  = SPREAD( zu, 2, statistic_regions+1 )
1737             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1738                dopr_initial_index(i) = 6
1739                hom(:,2,6,:)          = SPREAD( zu, 2, statistic_regions+1 )
1740                data_output_pr(i)     = data_output_pr(i)(2:)
1741             ENDIF
1742
1743          CASE ( 'w' )
1744             dopr_index(i) = 3
1745             dopr_unit(i)  = 'm/s'
1746             hom(:,2,3,:)  = SPREAD( zw, 2, statistic_regions+1 )
1747
1748          CASE ( 'theta', '#theta' )
1749             IF ( .NOT. bulk_cloud_model ) THEN
1750                dopr_index(i) = 4
1751                dopr_unit(i)  = 'K'
1752                hom(:,2,4,:)  = SPREAD( zu, 2, statistic_regions+1 )
1753                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1754                   dopr_initial_index(i) = 7
1755                   hom(:,2,7,:)          = SPREAD( zu, 2, statistic_regions+1 )
1756                   hom(nzb,2,7,:)        = 0.0_wp    ! because zu(nzb) is negative
1757                   data_output_pr(i)     = data_output_pr(i)(2:)
1758                ENDIF
1759             ELSE
1760                dopr_index(i) = 43
1761                dopr_unit(i)  = 'K'
1762                hom(:,2,43,:)  = SPREAD( zu, 2, statistic_regions+1 )
1763                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1764                   dopr_initial_index(i) = 28
1765                   hom(:,2,28,:)         = SPREAD( zu, 2, statistic_regions+1 )
1766                   hom(nzb,2,28,:)       = 0.0_wp    ! because zu(nzb) is negative
1767                   data_output_pr(i)     = data_output_pr(i)(2:)
1768                ENDIF
1769             ENDIF
1770
1771          CASE ( 'e', '#e' )
1772             dopr_index(i)  = 8
1773             dopr_unit(i)   = 'm2/s2'
1774             hom(:,2,8,:)   = SPREAD( zu, 2, statistic_regions+1 )
1775             hom(nzb,2,8,:) = 0.0_wp
1776             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1777                dopr_initial_index(i) = 8
1778                hom(:,2,8,:)          = SPREAD( zu, 2, statistic_regions+1 )
1779                data_output_pr(i)     = data_output_pr(i)(2:)
1780             ENDIF
1781
1782          CASE ( 'km', '#km' )
1783             dopr_index(i)  = 9
1784             dopr_unit(i)   = 'm2/s'
1785             hom(:,2,9,:)   = SPREAD( zu, 2, statistic_regions+1 )
1786             hom(nzb,2,9,:) = 0.0_wp
1787             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1788                dopr_initial_index(i) = 23
1789                hom(:,2,23,:)         = hom(:,2,9,:)
1790                data_output_pr(i)     = data_output_pr(i)(2:)
1791             ENDIF
1792
1793          CASE ( 'kh', '#kh' )
1794             dopr_index(i)   = 10
1795             dopr_unit(i)    = 'm2/s'
1796             hom(:,2,10,:)   = SPREAD( zu, 2, statistic_regions+1 )
1797             hom(nzb,2,10,:) = 0.0_wp
1798             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1799                dopr_initial_index(i) = 24
1800                hom(:,2,24,:)         = hom(:,2,10,:)
1801                data_output_pr(i)     = data_output_pr(i)(2:)
1802             ENDIF
1803
1804          CASE ( 'l', '#l' )
1805             dopr_index(i)   = 11
1806             dopr_unit(i)    = 'm'
1807             hom(:,2,11,:)   = SPREAD( zu, 2, statistic_regions+1 )
1808             hom(nzb,2,11,:) = 0.0_wp
1809             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1810                dopr_initial_index(i) = 25
1811                hom(:,2,25,:)         = hom(:,2,11,:)
1812                data_output_pr(i)     = data_output_pr(i)(2:)
1813             ENDIF
1814
1815          CASE ( 'w"u"' )
1816             dopr_index(i) = 12
1817             dopr_unit(i)  = TRIM ( momentumflux_output_unit )
1818             hom(:,2,12,:) = SPREAD( zw, 2, statistic_regions+1 )
1819             IF ( constant_flux_layer )  hom(nzb,2,12,:) = zu(1)
1820
1821          CASE ( 'w*u*' )
1822             dopr_index(i) = 13
1823             dopr_unit(i)  = TRIM ( momentumflux_output_unit )
1824             hom(:,2,13,:) = SPREAD( zw, 2, statistic_regions+1 )
1825
1826          CASE ( 'w"v"' )
1827             dopr_index(i) = 14
1828             dopr_unit(i)  = TRIM ( momentumflux_output_unit )
1829             hom(:,2,14,:) = SPREAD( zw, 2, statistic_regions+1 )
1830             IF ( constant_flux_layer )  hom(nzb,2,14,:) = zu(1)
1831
1832          CASE ( 'w*v*' )
1833             dopr_index(i) = 15
1834             dopr_unit(i)  = TRIM ( momentumflux_output_unit )
1835             hom(:,2,15,:) = SPREAD( zw, 2, statistic_regions+1 )
1836
1837          CASE ( 'w"theta"' )
1838             dopr_index(i) = 16
1839             dopr_unit(i)  = TRIM ( heatflux_output_unit )
1840             hom(:,2,16,:) = SPREAD( zw, 2, statistic_regions+1 )
1841
1842          CASE ( 'w*theta*' )
1843             dopr_index(i) = 17
1844             dopr_unit(i)  = TRIM ( heatflux_output_unit )
1845             hom(:,2,17,:) = SPREAD( zw, 2, statistic_regions+1 )
1846
1847          CASE ( 'wtheta' )
1848             dopr_index(i) = 18
1849             dopr_unit(i)  = TRIM ( heatflux_output_unit )
1850             hom(:,2,18,:) = SPREAD( zw, 2, statistic_regions+1 )
1851
1852          CASE ( 'wu' )
1853             dopr_index(i) = 19
1854             dopr_unit(i)  = TRIM ( momentumflux_output_unit )
1855             hom(:,2,19,:) = SPREAD( zw, 2, statistic_regions+1 )
1856             IF ( constant_flux_layer )  hom(nzb,2,19,:) = zu(1)
1857
1858          CASE ( 'wv' )
1859             dopr_index(i) = 20
1860             dopr_unit(i)  = TRIM ( momentumflux_output_unit )
1861             hom(:,2,20,:) = SPREAD( zw, 2, statistic_regions+1 )
1862             IF ( constant_flux_layer )  hom(nzb,2,20,:) = zu(1)
1863
1864          CASE ( 'w*theta*BC' )
1865             dopr_index(i) = 21
1866             dopr_unit(i)  = TRIM ( heatflux_output_unit )
1867             hom(:,2,21,:) = SPREAD( zw, 2, statistic_regions+1 )
1868
1869          CASE ( 'wthetaBC' )
1870             dopr_index(i) = 22
1871             dopr_unit(i)  = TRIM ( heatflux_output_unit )
1872             hom(:,2,22,:) = SPREAD( zw, 2, statistic_regions+1 )
1873
1874          CASE ( 'u*2' )
1875             dopr_index(i) = 30
1876             dopr_unit(i)  = 'm2/s2'
1877             hom(:,2,30,:) = SPREAD( zu, 2, statistic_regions+1 )
1878
1879          CASE ( 'v*2' )
1880             dopr_index(i) = 31
1881             dopr_unit(i)  = 'm2/s2'
1882             hom(:,2,31,:) = SPREAD( zu, 2, statistic_regions+1 )
1883
1884          CASE ( 'w*2' )
1885             dopr_index(i) = 32
1886             dopr_unit(i)  = 'm2/s2'
1887             hom(:,2,32,:) = SPREAD( zw, 2, statistic_regions+1 )
1888
1889          CASE ( 'theta*2' )
1890             dopr_index(i) = 33
1891             dopr_unit(i)  = 'K2'
1892             hom(:,2,33,:) = SPREAD( zu, 2, statistic_regions+1 )
1893
1894          CASE ( 'e*' )
1895             dopr_index(i) = 34
1896             dopr_unit(i)  = 'm2/s2'
1897             hom(:,2,34,:) = SPREAD( zu, 2, statistic_regions+1 )
1898
1899          CASE ( 'w*2theta*' )
1900             dopr_index(i) = 35
1901             dopr_unit(i)  = 'K m2/s2'
1902             hom(:,2,35,:) = SPREAD( zw, 2, statistic_regions+1 )
1903
1904          CASE ( 'w*theta*2' )
1905             dopr_index(i) = 36
1906             dopr_unit(i)  = 'K2 m/s'
1907             hom(:,2,36,:) = SPREAD( zw, 2, statistic_regions+1 )
1908
1909          CASE ( 'w*e*' )
1910             dopr_index(i) = 37
1911             dopr_unit(i)  = 'm3/s3'
1912             hom(:,2,37,:) = SPREAD( zw, 2, statistic_regions+1 )
1913
1914          CASE ( 'w*3' )
1915             dopr_index(i) = 38
1916             dopr_unit(i)  = 'm3/s3'
1917             hom(:,2,38,:) = SPREAD( zw, 2, statistic_regions+1 )
1918
1919          CASE ( 'Sw' )
1920             dopr_index(i) = 39
1921             dopr_unit(i)  = 'none'
1922             hom(:,2,39,:) = SPREAD( zw, 2, statistic_regions+1 )
1923
1924          CASE ( 'p' )
1925             dopr_index(i) = 40
1926             dopr_unit(i)  = 'Pa'
1927             hom(:,2,40,:) = SPREAD( zu, 2, statistic_regions+1 )
1928
1929          CASE ( 'q', '#q' )
1930             IF ( .NOT. humidity )  THEN
1931                message_string = 'data_output_pr = ' //                        &
1932                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
1933                                 'lemented for humidity = .FALSE.'
1934                CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 )
1935             ELSE
1936                dopr_index(i) = 41
1937                dopr_unit(i)  = 'kg/kg'
1938                hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 )
1939                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1940                   dopr_initial_index(i) = 26
1941                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1942                   hom(nzb,2,26,:)       = 0.0_wp    ! because zu(nzb) is negative
1943                   data_output_pr(i)     = data_output_pr(i)(2:)
1944                ENDIF
1945             ENDIF
1946
1947          CASE ( 's', '#s' )
1948             IF ( .NOT. passive_scalar )  THEN
1949                message_string = 'data_output_pr = ' //                        &
1950                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
1951                                 'lemented for passive_scalar = .FALSE.'
1952                CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 )
1953             ELSE
1954                dopr_index(i) = 115
1955                dopr_unit(i)  = 'kg/m3'
1956                hom(:,2,115,:) = SPREAD( zu, 2, statistic_regions+1 )
1957                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1958                   dopr_initial_index(i) = 121
1959                   hom(:,2,121,:)        = SPREAD( zu, 2, statistic_regions+1 )
1960                   hom(nzb,2,121,:)      = 0.0_wp    ! because zu(nzb) is negative
1961                   data_output_pr(i)     = data_output_pr(i)(2:)
1962                ENDIF
1963             ENDIF
1964
1965          CASE ( 'qv', '#qv' )
1966             IF ( .NOT. bulk_cloud_model ) THEN
1967                dopr_index(i) = 41
1968                dopr_unit(i)  = 'kg/kg'
1969                hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 )
1970                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1971                   dopr_initial_index(i) = 26
1972                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1973                   hom(nzb,2,26,:)       = 0.0_wp    ! because zu(nzb) is negative
1974                   data_output_pr(i)     = data_output_pr(i)(2:)
1975                ENDIF
1976             ELSE
1977                dopr_index(i) = 42
1978                dopr_unit(i)  = 'kg/kg'
1979                hom(:,2,42,:) = SPREAD( zu, 2, statistic_regions+1 )
1980                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1981                   dopr_initial_index(i) = 27
1982                   hom(:,2,27,:)         = SPREAD( zu, 2, statistic_regions+1 )
1983                   hom(nzb,2,27,:)       = 0.0_wp   ! because zu(nzb) is negative
1984                   data_output_pr(i)     = data_output_pr(i)(2:)
1985                ENDIF
1986             ENDIF
1987
1988          CASE ( 'thetal', '#thetal' )
1989             IF ( .NOT. bulk_cloud_model ) THEN
1990                message_string = 'data_output_pr = ' //                        &
1991                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
1992                                 'lemented for bulk_cloud_model = .FALSE.'
1993                CALL message( 'check_parameters', 'PA0094', 1, 2, 0, 6, 0 )
1994             ELSE
1995                dopr_index(i) = 4
1996                dopr_unit(i)  = 'K'
1997                hom(:,2,4,:)  = SPREAD( zu, 2, statistic_regions+1 )
1998                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1999                   dopr_initial_index(i) = 7
2000                   hom(:,2,7,:)          = SPREAD( zu, 2, statistic_regions+1 )
2001                   hom(nzb,2,7,:)        = 0.0_wp    ! because zu(nzb) is negative
2002                   data_output_pr(i)     = data_output_pr(i)(2:)
2003                ENDIF
2004             ENDIF
2005
2006          CASE ( 'thetav', '#thetav' )
2007             dopr_index(i) = 44
2008             dopr_unit(i)  = 'K'
2009             hom(:,2,44,:) = SPREAD( zu, 2, statistic_regions+1 )
2010             IF ( data_output_pr(i)(1:1) == '#' )  THEN
2011                dopr_initial_index(i) = 29
2012                hom(:,2,29,:)         = SPREAD( zu, 2, statistic_regions+1 )
2013                hom(nzb,2,29,:)       = 0.0_wp    ! because zu(nzb) is negative
2014                data_output_pr(i)     = data_output_pr(i)(2:)
2015             ENDIF
2016
2017          CASE ( 'w"thetav"' )
2018             dopr_index(i) = 45
2019             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2020             hom(:,2,45,:) = SPREAD( zw, 2, statistic_regions+1 )
2021
2022          CASE ( 'w*thetav*' )
2023             dopr_index(i) = 46
2024             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2025             hom(:,2,46,:) = SPREAD( zw, 2, statistic_regions+1 )
2026
2027          CASE ( 'wthetav' )
2028             dopr_index(i) = 47
2029             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2030             hom(:,2,47,:) = SPREAD( zw, 2, statistic_regions+1 )
2031
2032          CASE ( 'w"q"' )
2033             IF (  .NOT.  humidity )  THEN
2034                message_string = 'data_output_pr = ' //                        &
2035                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2036                                 'lemented for humidity = .FALSE.'
2037                CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 )
2038             ELSE
2039                dopr_index(i) = 48
2040                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2041                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
2042             ENDIF
2043
2044          CASE ( 'w*q*' )
2045             IF (  .NOT.  humidity )  THEN
2046                message_string = 'data_output_pr = ' //                        &
2047                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2048                                 'lemented for humidity = .FALSE.'
2049                CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 )
2050             ELSE
2051                dopr_index(i) = 49
2052                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2053                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
2054             ENDIF
2055
2056          CASE ( 'wq' )
2057             IF (  .NOT.  humidity )  THEN
2058                message_string = 'data_output_pr = ' //                        &
2059                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2060                                 'lemented for humidity = .FALSE.'
2061                CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 )
2062             ELSE
2063                dopr_index(i) = 50
2064                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2065                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
2066             ENDIF
2067
2068          CASE ( 'w"s"' )
2069             IF (  .NOT.  passive_scalar )  THEN
2070                message_string = 'data_output_pr = ' //                        &
2071                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2072                                 'lemented for passive_scalar = .FALSE.'
2073                CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 )
2074             ELSE
2075                dopr_index(i) = 117
2076                dopr_unit(i)  = 'kg/m3 m/s'
2077                hom(:,2,117,:) = SPREAD( zw, 2, statistic_regions+1 )
2078             ENDIF
2079
2080          CASE ( 'w*s*' )
2081             IF (  .NOT.  passive_scalar )  THEN
2082                message_string = 'data_output_pr = ' //                        &
2083                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2084                                 'lemented for passive_scalar = .FALSE.'
2085                CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 )
2086             ELSE
2087                dopr_index(i) = 114
2088                dopr_unit(i)  = 'kg/m3 m/s'
2089                hom(:,2,114,:) = SPREAD( zw, 2, statistic_regions+1 )
2090             ENDIF
2091
2092          CASE ( 'ws' )
2093             IF (  .NOT.  passive_scalar )  THEN
2094                message_string = 'data_output_pr = ' //                        &
2095                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2096                                 'lemented for passive_scalar = .FALSE.'
2097                CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 )
2098             ELSE
2099                dopr_index(i) = 118
2100                dopr_unit(i)  = 'kg/m3 m/s'
2101                hom(:,2,118,:) = SPREAD( zw, 2, statistic_regions+1 )
2102             ENDIF
2103
2104          CASE ( 'w"qv"' )
2105             IF ( humidity  .AND.  .NOT.  bulk_cloud_model )  THEN
2106                dopr_index(i) = 48
2107                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2108                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
2109             ELSEIF ( humidity  .AND.  bulk_cloud_model )  THEN
2110                dopr_index(i) = 51
2111                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2112                hom(:,2,51,:) = SPREAD( zw, 2, statistic_regions+1 )
2113             ELSE
2114                message_string = 'data_output_pr = ' //                        &
2115                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2116                                 'lemented for bulk_cloud_model = .FALSE. ' // &
2117                                 'and humidity = .FALSE.'
2118                CALL message( 'check_parameters', 'PA0095', 1, 2, 0, 6, 0 )
2119             ENDIF
2120
2121          CASE ( 'w*qv*' )
2122             IF ( humidity  .AND.  .NOT. bulk_cloud_model )  THEN
2123                dopr_index(i) = 49
2124                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2125                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
2126             ELSEIF( humidity .AND. bulk_cloud_model ) THEN
2127                dopr_index(i) = 52
2128                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2129                hom(:,2,52,:) = SPREAD( zw, 2, statistic_regions+1 )
2130             ELSE
2131                message_string = 'data_output_pr = ' //                        &
2132                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2133                                 'lemented for bulk_cloud_model = .FALSE. ' // &
2134                                 'and humidity = .FALSE.'
2135                CALL message( 'check_parameters', 'PA0095', 1, 2, 0, 6, 0 )
2136             ENDIF
2137
2138          CASE ( 'wqv' )
2139             IF ( humidity  .AND.  .NOT.  bulk_cloud_model )  THEN
2140                dopr_index(i) = 50
2141                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2142                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
2143             ELSEIF ( humidity  .AND.  bulk_cloud_model )  THEN
2144                dopr_index(i) = 53
2145                dopr_unit(i)  = TRIM ( waterflux_output_unit )
2146                hom(:,2,53,:) = SPREAD( zw, 2, statistic_regions+1 )
2147             ELSE
2148                message_string = 'data_output_pr = ' //                        &
2149                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2150                                 'lemented for bulk_cloud_model = .FALSE. ' // &
2151                                 'and humidity = .FALSE.'
2152                CALL message( 'check_parameters', 'PA0095', 1, 2, 0, 6, 0 )
2153             ENDIF
2154
2155          CASE ( 'ql' )
2156             IF (  .NOT.  bulk_cloud_model  .AND.  .NOT.  cloud_droplets )  THEN
2157                message_string = 'data_output_pr = ' //                        &
2158                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2159                                 'lemented for bulk_cloud_model = .FALSE. ' // &
2160                                 'and cloud_droplets = .FALSE.'
2161                CALL message( 'check_parameters', 'PA0096', 1, 2, 0, 6, 0 )
2162             ELSE
2163                dopr_index(i) = 54
2164                dopr_unit(i)  = 'kg/kg'
2165                hom(:,2,54,:)  = SPREAD( zu, 2, statistic_regions+1 )
2166             ENDIF
2167
2168          CASE ( 'w*u*u*:dz' )
2169             dopr_index(i) = 55
2170             dopr_unit(i)  = 'm2/s3'
2171             hom(:,2,55,:) = SPREAD( zu, 2, statistic_regions+1 )
2172
2173          CASE ( 'w*p*:dz' )
2174             dopr_index(i) = 56
2175             dopr_unit(i)  = 'm2/s3'
2176             hom(:,2,56,:) = SPREAD( zw, 2, statistic_regions+1 )
2177
2178          CASE ( 'w"e:dz' )
2179             dopr_index(i) = 57
2180             dopr_unit(i)  = 'm2/s3'
2181             hom(:,2,57,:) = SPREAD( zu, 2, statistic_regions+1 )
2182
2183          CASE ( 'u"theta"' )
2184             dopr_index(i) = 58
2185             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2186             hom(:,2,58,:) = SPREAD( zu, 2, statistic_regions+1 )
2187
2188          CASE ( 'u*theta*' )
2189             dopr_index(i) = 59
2190             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2191             hom(:,2,59,:) = SPREAD( zu, 2, statistic_regions+1 )
2192
2193          CASE ( 'utheta_t' )
2194             dopr_index(i) = 60
2195             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2196             hom(:,2,60,:) = SPREAD( zu, 2, statistic_regions+1 )
2197
2198          CASE ( 'v"theta"' )
2199             dopr_index(i) = 61
2200             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2201             hom(:,2,61,:) = SPREAD( zu, 2, statistic_regions+1 )
2202
2203          CASE ( 'v*theta*' )
2204             dopr_index(i) = 62
2205             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2206             hom(:,2,62,:) = SPREAD( zu, 2, statistic_regions+1 )
2207
2208          CASE ( 'vtheta_t' )
2209             dopr_index(i) = 63
2210             dopr_unit(i)  = TRIM ( heatflux_output_unit )
2211             hom(:,2,63,:) = SPREAD( zu, 2, statistic_regions+1 )
2212
2213          CASE ( 'w*p*' )
2214             dopr_index(i) = 68
2215             dopr_unit(i)  = 'm3/s3'
2216             hom(:,2,68,:) = SPREAD( zu, 2, statistic_regions+1 )
2217
2218          CASE ( 'w"e' )
2219             dopr_index(i) = 69
2220             dopr_unit(i)  = 'm3/s3'
2221             hom(:,2,69,:) = SPREAD( zu, 2, statistic_regions+1 )
2222
2223          CASE ( 'q*2' )
2224             IF (  .NOT.  humidity )  THEN
2225                message_string = 'data_output_pr = ' //                        &
2226                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2227                                 'lemented for humidity = .FALSE.'
2228                CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 )
2229             ELSE
2230                dopr_index(i) = 70
2231                dopr_unit(i)  = 'kg2/kg2'
2232                hom(:,2,70,:) = SPREAD( zu, 2, statistic_regions+1 )
2233             ENDIF
2234
2235          CASE ( 'hyp' )
2236             dopr_index(i) = 72
2237             dopr_unit(i)  = 'hPa'
2238             hom(:,2,72,:) = SPREAD( zu, 2, statistic_regions+1 )
2239
2240          CASE ( 'rho' )
2241             dopr_index(i)  = 119
2242             dopr_unit(i)   = 'kg/m3'
2243             hom(:,2,119,:) = SPREAD( zu, 2, statistic_regions+1 )
2244
2245          CASE ( 'rho_zw' )
2246             dopr_index(i)  = 120
2247             dopr_unit(i)   = 'kg/m3'
2248             hom(:,2,120,:) = SPREAD( zw, 2, statistic_regions+1 )
2249
2250          CASE ( 'ug' )
2251             dopr_index(i) = 78
2252             dopr_unit(i)  = 'm/s'
2253             hom(:,2,78,:) = SPREAD( zu, 2, statistic_regions+1 )
2254
2255          CASE ( 'vg' )
2256             dopr_index(i) = 79
2257             dopr_unit(i)  = 'm/s'
2258             hom(:,2,79,:) = SPREAD( zu, 2, statistic_regions+1 )
2259
2260          CASE ( 'w_subs' )
2261             IF (  .NOT.  large_scale_subsidence )  THEN
2262                message_string = 'data_output_pr = ' //                        &
2263                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2264                                 'lemented for large_scale_subsidence = .FALSE.'
2265                CALL message( 'check_parameters', 'PA0382', 1, 2, 0, 6, 0 )
2266             ELSE
2267                dopr_index(i) = 80
2268                dopr_unit(i)  = 'm/s'
2269                hom(:,2,80,:) = SPREAD( zu, 2, statistic_regions+1 )
2270             ENDIF
2271
2272          CASE ( 's*2' )
2273             IF (  .NOT.  passive_scalar )  THEN
2274                message_string = 'data_output_pr = ' //                        &
2275                                 TRIM( data_output_pr(i) ) // ' is not imp' // &
2276                                 'lemented for passive_scalar = .FALSE.'
2277                CALL message( 'check_parameters', 'PA0185', 1, 2, 0, 6, 0 )
2278             ELSE
2279                dopr_index(i) = 116
2280                dopr_unit(i)  = 'kg2/m6'
2281                hom(:,2,116,:) = SPREAD( zu, 2, statistic_regions+1 )
2282             ENDIF
2283
2284           CASE ( 'eta' )
2285              dopr_index(i) = 121
2286              dopr_unit(i)  = 'mm'
2287              hom(:,2,121,:) = SPREAD( zu, 2, statistic_regions+1 )
2288
2289              kolmogorov_length_scale = .TRUE.
2290
2291          CASE DEFAULT
2292             unit = 'illegal'
2293!
2294!--          Check for other modules
2295             CALL module_interface_check_data_output_pr( data_output_pr(i), i, &
2296                                                         unit, dopr_unit(i) )
2297
2298!
2299!--          No valid quantity found
2300             IF ( unit == 'illegal' )  THEN
2301                IF ( data_output_pr_user(1) /= ' ' )  THEN
2302                   message_string = 'illegal value for data_output_pr or ' //  &
2303                                    'data_output_pr_user = "' //               &
2304                                    TRIM( data_output_pr(i) ) // '"'
2305                   CALL message( 'check_parameters', 'PA0097', 1, 2, 0, 6, 0 )
2306                ELSE
2307                   message_string = 'illegal value for data_output_pr = "' //  &
2308                                    TRIM( data_output_pr(i) ) // '"'
2309                   CALL message( 'check_parameters', 'PA0098', 1, 2, 0, 6, 0 )
2310                ENDIF
2311             ENDIF
2312
2313       END SELECT
2314
2315    ENDDO
2316
2317
2318!
2319!-- Append user-defined data output variables to the standard data output
2320    IF ( data_output_user(1) /= ' ' )  THEN
2321       i = 1
2322       DO  WHILE ( data_output(i) /= ' '  .AND.  i <= 500 )
2323          i = i + 1
2324       ENDDO
2325       j = 1
2326       DO  WHILE ( data_output_user(j) /= ' '  .AND.  j <= 500 )
2327          IF ( i > 500 )  THEN
2328             message_string = 'number of output quantitities given by data' // &
2329                '_output and data_output_user exceeds the limit of 500'
2330             CALL message( 'check_parameters', 'PA0102', 1, 2, 0, 6, 0 )
2331          ENDIF
2332          data_output(i) = data_output_user(j)
2333          i = i + 1
2334          j = j + 1
2335       ENDDO
2336    ENDIF
2337
2338!
2339!-- Check and set steering parameters for 2d/3d data output and averaging
2340    i   = 1
2341    DO  WHILE ( data_output(i) /= ' '  .AND.  i <= 500 )
2342!
2343!--    Check for data averaging
2344       ilen = LEN_TRIM( data_output(i) )
2345       j = 0                                                 ! no data averaging
2346       IF ( ilen > 3 )  THEN
2347          IF ( data_output(i)(ilen-2:ilen) == '_av' )  THEN
2348             j = 1                                           ! data averaging
2349             data_output(i) = data_output(i)(1:ilen-3)
2350          ENDIF
2351       ENDIF
2352!
2353!--    Check for cross section or volume data
2354       ilen = LEN_TRIM( data_output(i) )
2355       k = 0                                                   ! 3d data
2356       var = data_output(i)(1:ilen)
2357       IF ( ilen > 3 )  THEN
2358          IF ( data_output(i)(ilen-2:ilen) == '_xy'  .OR.                      &
2359               data_output(i)(ilen-2:ilen) == '_xz'  .OR.                      &
2360               data_output(i)(ilen-2:ilen) == '_yz' )  THEN
2361             k = 1                                             ! 2d data
2362             var = data_output(i)(1:ilen-3)
2363          ENDIF
2364       ENDIF
2365
2366!
2367!--    Check for allowed value and set units
2368       SELECT CASE ( TRIM( var ) )
2369
2370          CASE ( 'e' )
2371             IF ( constant_diffusion )  THEN
2372                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2373                                 'res constant_diffusion = .FALSE.'
2374                CALL message( 'check_parameters', 'PA0103', 1, 2, 0, 6, 0 )
2375             ENDIF
2376             unit = 'm2/s2'
2377
2378          CASE ( 'thetal' )
2379             IF (  .NOT.  bulk_cloud_model )  THEN
2380                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2381                         'res bulk_cloud_model = .TRUE.'
2382                CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 )
2383             ENDIF
2384             unit = 'K'
2385
2386          CASE ( 'pc', 'pr' )
2387             IF (  .NOT.  particle_advection )  THEN
2388                message_string = 'output of "' // TRIM( var ) // '" requir' // &
2389                   'es a "particle_parameters"-NAMELIST in the parameter ' //  &
2390                   'file (PARIN)'
2391                CALL message( 'check_parameters', 'PA0104', 1, 2, 0, 6, 0 )
2392             ENDIF
2393             IF ( TRIM( var ) == 'pc' )  unit = 'number'
2394             IF ( TRIM( var ) == 'pr' )  unit = 'm'
2395
2396          CASE ( 'q', 'thetav' )
2397             IF (  .NOT.  humidity )  THEN
2398                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2399                                 'res humidity = .TRUE.'
2400                CALL message( 'check_parameters', 'PA0105', 1, 2, 0, 6, 0 )
2401             ENDIF
2402             IF ( TRIM( var ) == 'q'   )  unit = 'kg/kg'
2403             IF ( TRIM( var ) == 'thetav' )  unit = 'K'
2404
2405          CASE ( 'ql' )
2406             IF ( .NOT.  ( bulk_cloud_model  .OR.  cloud_droplets ) )  THEN
2407                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2408                      'res bulk_cloud_model = .TRUE. or cloud_droplets = .TRUE.'
2409                CALL message( 'check_parameters', 'PA0106', 1, 2, 0, 6, 0 )
2410             ENDIF
2411             unit = 'kg/kg'
2412
2413          CASE ( 'ql_c', 'ql_v', 'ql_vp' )
2414             IF (  .NOT.  cloud_droplets )  THEN
2415                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2416                                 'res cloud_droplets = .TRUE.'
2417                CALL message( 'check_parameters', 'PA0107', 1, 2, 0, 6, 0 )
2418             ENDIF
2419             IF ( TRIM( var ) == 'ql_c'  )  unit = 'kg/kg'
2420             IF ( TRIM( var ) == 'ql_v'  )  unit = 'm3'
2421             IF ( TRIM( var ) == 'ql_vp' )  unit = 'none'
2422
2423          CASE ( 'qv' )
2424             IF (  .NOT.  bulk_cloud_model )  THEN
2425                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2426                                 'res bulk_cloud_model = .TRUE.'
2427                CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 )
2428             ENDIF
2429             unit = 'kg/kg'
2430
2431          CASE ( 's' )
2432             IF (  .NOT.  passive_scalar )  THEN
2433                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2434                                 'res passive_scalar = .TRUE.'
2435                CALL message( 'check_parameters', 'PA0110', 1, 2, 0, 6, 0 )
2436             ENDIF
2437             unit = 'kg/m3'
2438
2439          CASE ( 'p', 'theta', 'u', 'v', 'w' )
2440             IF ( TRIM( var ) == 'p'  )  unit = 'Pa'
2441             IF ( TRIM( var ) == 'theta' )  unit = 'K'
2442             IF ( TRIM( var ) == 'u'  )  unit = 'm/s'
2443             IF ( TRIM( var ) == 'v'  )  unit = 'm/s'
2444             IF ( TRIM( var ) == 'w'  )  unit = 'm/s'
2445             CONTINUE
2446
2447          CASE ( 'ghf*', 'lwp*', 'ol*', 'qsws*', 'r_a*',                       &
2448                 'shf*', 'ssws*', 't*', 'tsurf*', 'us*',                       &
2449                 'z0*', 'z0h*', 'z0q*' )
2450             IF ( k == 0  .OR.  data_output(i)(ilen-2:ilen) /= '_xy' )  THEN
2451                message_string = 'illegal value for data_output: "' //         &
2452                                 TRIM( var ) // '" & only 2d-horizontal ' //   &
2453                                 'cross sections are allowed for this value'
2454                CALL message( 'check_parameters', 'PA0111', 1, 2, 0, 6, 0 )
2455             ENDIF
2456
2457             IF ( TRIM( var ) == 'lwp*'  .AND.  .NOT. bulk_cloud_model )  THEN
2458                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2459                                 'res bulk_cloud_model = .TRUE.'
2460                CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 )
2461             ENDIF
2462             IF ( TRIM( var ) == 'qsws*'  .AND.  .NOT.  humidity )  THEN
2463                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2464                                 'res humidity = .TRUE.'
2465                CALL message( 'check_parameters', 'PA0322', 1, 2, 0, 6, 0 )
2466             ENDIF
2467
2468             IF ( TRIM( var ) == 'ghf*'  .AND.  .NOT.  land_surface )  THEN
2469                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2470                                 'res land_surface = .TRUE.'
2471                CALL message( 'check_parameters', 'PA0404', 1, 2, 0, 6, 0 )
2472             ENDIF
2473
2474             IF ( ( TRIM( var ) == 'r_a*' .OR.  TRIM( var ) == 'ghf*' )        &
2475                 .AND.  .NOT.  land_surface  .AND.  .NOT.  urban_surface )     &
2476             THEN
2477                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2478                                 'res land_surface = .TRUE. or ' //            &
2479                                 'urban_surface = .TRUE.'
2480                CALL message( 'check_parameters', 'PA0404', 1, 2, 0, 6, 0 )
2481             ENDIF
2482
2483             IF ( TRIM( var ) == 'ssws*'  .AND.  .NOT.  passive_scalar )  THEN
2484                message_string = 'output of "' // TRIM( var ) // '" requi' //  &
2485                                 'res passive_scalar = .TRUE.'
2486                CALL message( 'check_parameters', 'PA0361', 1, 2, 0, 6, 0 )
2487             ENDIF
2488
2489             IF ( TRIM( var ) == 'ghf*'   )  unit = 'W/m2'
2490             IF ( TRIM( var ) == 'lwp*'   )  unit = 'kg/m2'
2491             IF ( TRIM( var ) == 'ol*'    )  unit = 'm'
2492             IF ( TRIM( var ) == 'qsws*'  )  unit = 'kgm/kgs'
2493             IF ( TRIM( var ) == 'r_a*'   )  unit = 's/m'
2494             IF ( TRIM( var ) == 'shf*'   )  unit = 'K*m/s'
2495             IF ( TRIM( var ) == 'ssws*'  )  unit = 'kg/m2*s'
2496             IF ( TRIM( var ) == 't*'     )  unit = 'K'
2497             IF ( TRIM( var ) == 'tsurf*' )  unit = 'K'
2498             IF ( TRIM( var ) == 'us*'    )  unit = 'm/s'
2499             IF ( TRIM( var ) == 'z0*'    )  unit = 'm'
2500             IF ( TRIM( var ) == 'z0h*'   )  unit = 'm'
2501!
2502!--          Output of surface latent and sensible heat flux will be in W/m2
2503!--          in case of natural- and urban-type surfaces, even if
2504!--          flux_output_mode is set to kinematic units.
2505             IF ( land_surface  .OR.  urban_surface )  THEN
2506                IF ( TRIM( var ) == 'shf*'   )  unit = 'W/m2'
2507                IF ( TRIM( var ) == 'qsws*'  )  unit = 'W/m2'
2508             ENDIF
2509
2510          CASE DEFAULT
2511!
2512!--          Check for other modules
2513             CALL module_interface_check_data_output( var, unit, i, j, ilen, k )
2514
2515             IF ( unit == 'illegal' )  THEN
2516                IF ( data_output_user(1) /= ' ' )  THEN
2517                   message_string = 'illegal value for data_output or ' //     &
2518                         'data_output_user = "' // TRIM( data_output(i) ) // '"'
2519                   CALL message( 'check_parameters', 'PA0114', 1, 2, 0, 6, 0 )
2520                ELSE
2521                   message_string = 'illegal value for data_output = "' //     &
2522                                    TRIM( data_output(i) ) // '"'
2523                   CALL message( 'check_parameters', 'PA0115', 1, 2, 0, 6, 0 )
2524                ENDIF
2525             ENDIF
2526
2527       END SELECT
2528!
2529!--    Set the internal steering parameters appropriately
2530       IF ( k == 0 )  THEN
2531          do3d_no(j)              = do3d_no(j) + 1
2532          do3d(j,do3d_no(j))      = data_output(i)
2533          do3d_unit(j,do3d_no(j)) = unit
2534       ELSE
2535          do2d_no(j)              = do2d_no(j) + 1
2536          do2d(j,do2d_no(j))      = data_output(i)
2537          do2d_unit(j,do2d_no(j)) = unit
2538          IF ( data_output(i)(ilen-2:ilen) == '_xy' )  THEN
2539             data_output_xy(j) = .TRUE.
2540          ENDIF
2541          IF ( data_output(i)(ilen-2:ilen) == '_xz' )  THEN
2542             data_output_xz(j) = .TRUE.
2543          ENDIF
2544          IF ( data_output(i)(ilen-2:ilen) == '_yz' )  THEN
2545             data_output_yz(j) = .TRUE.
2546          ENDIF
2547       ENDIF
2548
2549       IF ( j == 1 )  THEN
2550!
2551!--       Check, if variable is already subject to averaging
2552          found = .FALSE.
2553          DO  k = 1, doav_n
2554             IF ( TRIM( doav(k) ) == TRIM( var ) )  found = .TRUE.
2555          ENDDO
2556
2557          IF ( .NOT. found )  THEN
2558             doav_n = doav_n + 1
2559             doav(doav_n) = var
2560          ENDIF
2561       ENDIF
2562
2563       i = i + 1
2564    ENDDO
2565
2566!
2567!-- Averaged 2d or 3d output requires that an averaging interval has been set
2568    IF ( doav_n > 0  .AND.  averaging_interval == 0.0_wp )  THEN
2569       WRITE( message_string, * )  'output of averaged quantity "',            &
2570                                   TRIM( doav(1) ), '_av" requires to set a ', &
2571                                   'non-zero averaging interval'
2572       CALL message( 'check_parameters', 'PA0323', 1, 2, 0, 6, 0 )
2573    ENDIF
2574
2575!
2576!-- Check sectional planes and store them in one shared array
2577    IF ( ANY( section_xy > nz + 1 ) )  THEN
2578       WRITE( message_string, * )  'section_xy must be <= nz + 1 = ', nz + 1
2579       CALL message( 'check_parameters', 'PA0319', 1, 2, 0, 6, 0 )
2580    ENDIF
2581    IF ( ANY( section_xz > ny + 1 ) )  THEN
2582       WRITE( message_string, * )  'section_xz must be <= ny + 1 = ', ny + 1
2583       CALL message( 'check_parameters', 'PA0320', 1, 2, 0, 6, 0 )
2584    ENDIF
2585    IF ( ANY( section_yz > nx + 1 ) )  THEN
2586       WRITE( message_string, * )  'section_yz must be <= nx + 1 = ', nx + 1
2587       CALL message( 'check_parameters', 'PA0321', 1, 2, 0, 6, 0 )
2588    ENDIF
2589    section(:,1) = section_xy
2590    section(:,2) = section_xz
2591    section(:,3) = section_yz
2592
2593    IF ( ANY( data_output_xy ) .AND. .NOT. ANY( section(:,1) /= -9999 ) )  THEN
2594       WRITE( message_string, * )  'section_xy not defined for requested '  // &
2595                                   'xy-cross section output.&At least one ' // &
2596                                   'cross section must be given.'
2597       CALL message( 'check_parameters', 'PA0681', 1, 2, 0, 6, 0 )
2598    ENDIF
2599    IF ( ANY( data_output_xz ) .AND. .NOT. ANY( section(:,2) /= -9999 ) )  THEN
2600       WRITE( message_string, * )  'section_xz not defined for requested '  // &
2601                                   'xz-cross section output.&At least one ' // &
2602                                   'cross section must be given.'
2603       CALL message( 'check_parameters', 'PA0681', 1, 2, 0, 6, 0 )
2604    ENDIF
2605    IF ( ANY( data_output_yz ) .AND. .NOT. ANY( section(:,3) /= -9999 ) )  THEN
2606       WRITE( message_string, * )  'section_yz not defined for requested '  // &
2607                                   'yz-cross section output.&At least one ' // &
2608                                   'cross section must be given.'
2609       CALL message( 'check_parameters', 'PA0681', 1, 2, 0, 6, 0 )
2610    ENDIF
2611!
2612!-- Upper plot limit for 3D arrays
2613    IF ( nz_do3d == -9999 )  nz_do3d = nzt + 1
2614
2615!
2616!-- Set output format string (used in header)
2617    SELECT CASE ( netcdf_data_format )
2618       CASE ( 1 )
2619          netcdf_data_format_string = 'netCDF classic'
2620       CASE ( 2 )
2621          netcdf_data_format_string = 'netCDF 64bit offset'
2622       CASE ( 3 )
2623          netcdf_data_format_string = 'netCDF4/HDF5'
2624       CASE ( 4 )
2625          netcdf_data_format_string = 'netCDF4/HDF5 classic'
2626       CASE ( 5 )
2627          netcdf_data_format_string = 'parallel netCDF4/HDF5'
2628       CASE ( 6 )
2629          netcdf_data_format_string = 'parallel netCDF4/HDF5 classic'
2630
2631    END SELECT
2632
2633!
2634!-- Check mask conditions
2635    DO mid = 1, max_masks
2636       IF ( data_output_masks(mid,1) /= ' '  .OR.                              &
2637            data_output_masks_user(mid,1) /= ' ' ) THEN
2638          masks = masks + 1
2639       ENDIF
2640    ENDDO
2641
2642    IF ( masks < 0  .OR.  masks > max_masks )  THEN
2643       WRITE( message_string, * )  'illegal value: masks must be >= 0 and ',   &
2644            '<= ', max_masks, ' (=max_masks)'
2645       CALL message( 'check_parameters', 'PA0325', 1, 2, 0, 6, 0 )
2646    ENDIF
2647    IF ( masks > 0 )  THEN
2648       mask_scale(1) = mask_scale_x
2649       mask_scale(2) = mask_scale_y
2650       mask_scale(3) = mask_scale_z
2651       IF ( ANY( mask_scale <= 0.0_wp ) )  THEN
2652          WRITE( message_string, * )                                           &
2653               'illegal value: mask_scale_x, mask_scale_y and mask_scale_z',   &
2654               'must be > 0.0'
2655          CALL message( 'check_parameters', 'PA0326', 1, 2, 0, 6, 0 )
2656       ENDIF
2657!
2658!--    Generate masks for masked data output
2659!--    Parallel netcdf output is not tested so far for masked data, hence
2660!--    netcdf_data_format is switched back to non-parallel output.
2661       netcdf_data_format_save = netcdf_data_format
2662       IF ( netcdf_data_format > 4 )  THEN
2663          IF ( netcdf_data_format == 5 ) netcdf_data_format = 3
2664          IF ( netcdf_data_format == 6 ) netcdf_data_format = 4
2665          message_string = 'netCDF file formats '//                            &
2666                           '5 (parallel netCDF 4) and ' //                     &
2667                           '6 (parallel netCDF 4 Classic model) '//            &
2668                           '& are currently not supported (not yet tested) ' //&
2669                           'for masked data. &Using respective non-parallel' //&
2670                           ' output for masked data.'
2671          CALL message( 'check_parameters', 'PA0383', 0, 0, 0, 6, 0 )
2672       ENDIF
2673       CALL init_masks
2674       netcdf_data_format = netcdf_data_format_save
2675    ENDIF
2676
2677!
2678!-- Check the NetCDF data format
2679    IF ( netcdf_data_format > 2 )  THEN
2680#if defined( __netcdf4 )
2681       CONTINUE
2682#else
2683       message_string = 'netCDF: netCDF4 format requested but no ' //          &
2684                        'cpp-directive __netcdf4 given & switch '  //          &
2685                        'back to 64-bit offset format'
2686       CALL message( 'check_parameters', 'PA0171', 0, 1, 0, 6, 0 )
2687       netcdf_data_format = 2
2688#endif
2689    ENDIF
2690    IF ( netcdf_data_format > 4 )  THEN
2691#if defined( __netcdf4 ) && defined( __netcdf4_parallel )
2692       CONTINUE
2693#else
2694       message_string = 'netCDF: netCDF4 parallel output requested but no ' // &
2695                        'cpp-directive __netcdf4_parallel given & switch '   //&
2696                        'back to netCDF4 non-parallel output'
2697       CALL message( 'check_parameters', 'PA0099', 0, 1, 0, 6, 0 )
2698       netcdf_data_format = netcdf_data_format - 2
2699#endif
2700    ENDIF
2701
2702!
2703!-- Calculate fixed number of output time levels for parallel netcdf output.
2704!-- The time dimension has to be defined as limited for parallel output,
2705!-- because otherwise the I/O performance drops significantly.
2706    IF ( netcdf_data_format > 4 )  THEN
2707
2708!
2709!--    Check if any of the follwoing data output interval is 0.0s, which is
2710!--    not allowed for parallel output.
2711       CALL check_dt_do( dt_do3d,           'dt_do3d'           )
2712       CALL check_dt_do( dt_do2d_xy,        'dt_do2d_xy'        )
2713       CALL check_dt_do( dt_do2d_xz,        'dt_do2d_xz'        )
2714       CALL check_dt_do( dt_do2d_yz,        'dt_do2d_yz'        )
2715       CALL check_dt_do( dt_data_output_av, 'dt_data_output_av' )
2716
2717!--    Set needed time levels (ntdim) to
2718!--    saved time levels + to be saved time levels.
2719       ntdim_3d(0) = do3d_time_count(0) + CEILING(                                    &
2720                     ( end_time - MAX(                                                &
2721                         MERGE(skip_time_do3d, skip_time_do3d + spinup_time,          &
2722                               data_output_during_spinup ),                           &
2723                         simulated_time_at_begin )                                    &
2724                     ) / dt_do3d )
2725       IF ( do3d_at_begin ) ntdim_3d(0) = ntdim_3d(0) + 1
2726
2727       ntdim_3d(1) = do3d_time_count(1) + CEILING(                                    &
2728                     ( end_time - MAX(                                                &
2729                         MERGE(   skip_time_data_output_av, skip_time_data_output_av  &
2730                                + spinup_time, data_output_during_spinup ),           &
2731                         simulated_time_at_begin )                                    &
2732                     ) / dt_data_output_av )
2733
2734       ntdim_2d_xy(0) = do2d_xy_time_count(0) + CEILING(                              &
2735                        ( end_time - MAX(                                             &
2736                           MERGE(skip_time_do2d_xy, skip_time_do2d_xy + spinup_time,  &
2737                                 data_output_during_spinup ),                         &
2738                           simulated_time_at_begin )                                  &
2739                        ) / dt_do2d_xy )
2740
2741       ntdim_2d_xz(0) = do2d_xz_time_count(0) + CEILING(                              &
2742                        ( end_time - MAX(                                             &
2743                         MERGE(skip_time_do2d_xz, skip_time_do2d_xz + spinup_time,    &
2744                               data_output_during_spinup ),                           &
2745                         simulated_time_at_begin )                                    &
2746                        ) / dt_do2d_xz )
2747
2748       ntdim_2d_yz(0) = do2d_yz_time_count(0) + CEILING(                              &
2749                        ( end_time - MAX(                                             &
2750                         MERGE(skip_time_do2d_yz, skip_time_do2d_yz + spinup_time,    &
2751                               data_output_during_spinup ),                           &
2752                         simulated_time_at_begin )                                    &
2753                        ) / dt_do2d_yz )
2754
2755       IF ( do2d_at_begin )  THEN
2756          ntdim_2d_xy(0) = ntdim_2d_xy(0) + 1
2757          ntdim_2d_xz(0) = ntdim_2d_xz(0) + 1
2758          ntdim_2d_yz(0) = ntdim_2d_yz(0) + 1
2759       ENDIF
2760!
2761!--    Please note, for averaged 2D data skip_time_data_output_av is the relavant
2762!--    output control parameter.
2763       ntdim_2d_xy(1) = do2d_xy_time_count(1) + CEILING(                              &
2764                     ( end_time - MAX( MERGE( skip_time_data_output_av,               &
2765                                              skip_time_data_output_av + spinup_time, &
2766                                              data_output_during_spinup ),            &
2767                                       simulated_time_at_begin )                      &
2768                     ) / dt_data_output_av )
2769
2770       ntdim_2d_xz(1) = do2d_xz_time_count(1) + CEILING(                              &
2771                     ( end_time - MAX( MERGE( skip_time_data_output_av,               &
2772                                              skip_time_data_output_av + spinup_time, &
2773                                              data_output_during_spinup ),            &
2774                                       simulated_time_at_begin )                      &
2775                     ) / dt_data_output_av )
2776
2777       ntdim_2d_yz(1) = do2d_yz_time_count(1) + CEILING(                              &
2778                     ( end_time - MAX( MERGE( skip_time_data_output_av,               &
2779                                              skip_time_data_output_av + spinup_time, &
2780                                              data_output_during_spinup ),            &
2781                                       simulated_time_at_begin )                      &
2782                     ) / dt_data_output_av )
2783
2784    ENDIF
2785
2786!
2787!-- Check, whether a constant diffusion coefficient shall be used
2788    IF ( km_constant /= -1.0_wp )  THEN
2789       IF ( km_constant < 0.0_wp )  THEN
2790          WRITE( message_string, * )  'km_constant = ', km_constant, ' < 0.0'
2791          CALL message( 'check_parameters', 'PA0121', 1, 2, 0, 6, 0 )
2792       ELSE
2793          IF ( prandtl_number < 0.0_wp )  THEN
2794             WRITE( message_string, * )  'prandtl_number = ', prandtl_number,  &
2795                                         ' < 0.0'
2796             CALL message( 'check_parameters', 'PA0122', 1, 2, 0, 6, 0 )
2797          ENDIF
2798          constant_diffusion = .TRUE.
2799
2800          IF ( constant_flux_layer )  THEN
2801             message_string = 'constant_flux_layer is not allowed with fixed ' &
2802                              // 'value of km'
2803             CALL message( 'check_parameters', 'PA0123', 1, 2, 0, 6, 0 )
2804          ENDIF
2805       ENDIF
2806    ENDIF
2807
2808!
2809!-- In case of non-cyclic lateral boundaries and a damping layer for the
2810!-- potential temperature, check the width of the damping layer
2811    IF ( bc_lr /= 'cyclic' ) THEN
2812       IF ( pt_damping_width < 0.0_wp  .OR.                                    &
2813            pt_damping_width > REAL( (nx+1) * dx ) )  THEN
2814          message_string = 'pt_damping_width out of range'
2815          CALL message( 'check_parameters', 'PA0124', 1, 2, 0, 6, 0 )
2816       ENDIF
2817    ENDIF
2818
2819    IF ( bc_ns /= 'cyclic' )  THEN
2820       IF ( pt_damping_width < 0.0_wp  .OR.                                    &
2821            pt_damping_width > REAL( (ny+1) * dy ) )  THEN
2822          message_string = 'pt_damping_width out of range'
2823          CALL message( 'check_parameters', 'PA0124', 1, 2, 0, 6, 0 )
2824       ENDIF
2825    ENDIF
2826
2827!
2828!-- Check value range for zeta = z/L
2829    IF ( zeta_min >= zeta_max )  THEN
2830       WRITE( message_string, * )  'zeta_min = ', zeta_min, ' must be less ',  &
2831                                   'than zeta_max = ', zeta_max
2832       CALL message( 'check_parameters', 'PA0125', 1, 2, 0, 6, 0 )
2833    ENDIF
2834
2835!
2836!-- Check random generator
2837    IF ( (random_generator /= 'system-specific'      .AND.                     &
2838          random_generator /= 'random-parallel'   )  .AND.                     &
2839          random_generator /= 'numerical-recipes' )  THEN
2840       message_string = 'unknown random generator: random_generator = "' //    &
2841                        TRIM( random_generator ) // '"'
2842       CALL message( 'check_parameters', 'PA0135', 1, 2, 0, 6, 0 )
2843    ENDIF
2844
2845!
2846!-- Determine upper and lower hight level indices for random perturbations
2847    IF ( disturbance_level_b == -9999999.9_wp )  THEN
2848       IF ( ocean_mode )  THEN
2849          disturbance_level_b     = zu((nzt*2)/3)
2850          disturbance_level_ind_b = ( nzt * 2 ) / 3
2851       ELSE
2852          disturbance_level_b     = zu(nzb+3)
2853          disturbance_level_ind_b = nzb + 3
2854       ENDIF
2855    ELSEIF ( disturbance_level_b < zu(3) )  THEN
2856       WRITE( message_string, * )  'disturbance_level_b = ',                   &
2857                           disturbance_level_b, ' must be >= ', zu(3), '(zu(3))'
2858       CALL message( 'check_parameters', 'PA0126', 1, 2, 0, 6, 0 )
2859    ELSEIF ( disturbance_level_b > zu(nzt-2) )  THEN
2860       WRITE( message_string, * )  'disturbance_level_b = ',                   &
2861                   disturbance_level_b, ' must be <= ', zu(nzt-2), '(zu(nzt-2))'
2862       CALL message( 'check_parameters', 'PA0127', 1, 2, 0, 6, 0 )
2863    ELSE
2864       DO  k = 3, nzt-2
2865          IF ( disturbance_level_b <= zu(k) )  THEN
2866             disturbance_level_ind_b = k
2867             EXIT
2868          ENDIF
2869       ENDDO
2870    ENDIF
2871
2872    IF ( disturbance_level_t == -9999999.9_wp )  THEN
2873       IF ( ocean_mode )  THEN
2874          disturbance_level_t     = zu(nzt-3)
2875          disturbance_level_ind_t = nzt - 3
2876       ELSE
2877          disturbance_level_t     = zu(nzt/3)
2878          disturbance_level_ind_t = nzt / 3
2879       ENDIF
2880    ELSEIF ( disturbance_level_t > zu(nzt-2) )  THEN
2881       WRITE( message_string, * )  'disturbance_level_t = ',                   &
2882                   disturbance_level_t, ' must be <= ', zu(nzt-2), '(zu(nzt-2))'
2883       CALL message( 'check_parameters', 'PA0128', 1, 2, 0, 6, 0 )
2884    ELSEIF ( disturbance_level_t < disturbance_level_b )  THEN
2885       WRITE( message_string, * )  'disturbance_level_t = ',                   &
2886                   disturbance_level_t, ' must be >= disturbance_level_b = ',  &
2887                   disturbance_level_b
2888       CALL message( 'check_parameters', 'PA0129', 1, 2, 0, 6, 0 )
2889    ELSE
2890       DO  k = 3, nzt-2
2891          IF ( disturbance_level_t <= zu(k) )  THEN
2892             disturbance_level_ind_t = k
2893             EXIT
2894          ENDIF
2895       ENDDO
2896    ENDIF
2897
2898!
2899!-- Check again whether the levels determined this way are ok.
2900!-- Error may occur at automatic determination and too few grid points in
2901!-- z-direction.
2902    IF ( disturbance_level_ind_t < disturbance_level_ind_b )  THEN
2903       WRITE( message_string, * )  'disturbance_level_ind_t = ',               &
2904                disturbance_level_ind_t, ' must be >= ',                       &
2905                'disturbance_level_ind_b = ', disturbance_level_ind_b
2906       CALL message( 'check_parameters', 'PA0130', 1, 2, 0, 6, 0 )
2907    ENDIF
2908
2909!
2910!-- Determine the horizontal index range for random perturbations.
2911!-- In case of non-cyclic horizontal boundaries, no perturbations are imposed
2912!-- near the inflow and the perturbation area is further limited to ...(1)
2913!-- after the initial phase of the flow.
2914
2915    IF ( bc_lr /= 'cyclic' )  THEN
2916       IF ( inflow_disturbance_begin == -1 )  THEN
2917          inflow_disturbance_begin = MIN( 10, nx/2 )
2918       ENDIF
2919       IF ( inflow_disturbance_begin < 0  .OR.  inflow_disturbance_begin > nx )&
2920       THEN
2921          message_string = 'inflow_disturbance_begin out of range'
2922          CALL message( 'check_parameters', 'PA0131', 1, 2, 0, 6, 0 )
2923       ENDIF
2924       IF ( inflow_disturbance_end == -1 )  THEN
2925          inflow_disturbance_end = MIN( 100, 3*nx/4 )
2926       ENDIF
2927       IF ( inflow_disturbance_end < 0  .OR.  inflow_disturbance_end > nx )    &
2928       THEN
2929          message_string = 'inflow_disturbance_end out of range'
2930          CALL message( 'check_parameters', 'PA0132', 1, 2, 0, 6, 0 )
2931       ENDIF
2932    ELSEIF ( bc_ns /= 'cyclic' )  THEN
2933       IF ( inflow_disturbance_begin == -1 )  THEN
2934          inflow_disturbance_begin = MIN( 10, ny/2 )
2935       ENDIF
2936       IF ( inflow_disturbance_begin < 0  .OR.  inflow_disturbance_begin > ny )&
2937       THEN
2938          message_string = 'inflow_disturbance_begin out of range'
2939          CALL message( 'check_parameters', 'PA0131', 1, 2, 0, 6, 0 )
2940       ENDIF
2941       IF ( inflow_disturbance_end == -1 )  THEN
2942          inflow_disturbance_end = MIN( 100, 3*ny/4 )
2943       ENDIF
2944       IF ( inflow_disturbance_end < 0  .OR.  inflow_disturbance_end > ny )    &
2945       THEN
2946          message_string = 'inflow_disturbance_end out of range'
2947          CALL message( 'check_parameters', 'PA0132', 1, 2, 0, 6, 0 )
2948       ENDIF
2949    ENDIF
2950
2951    IF ( random_generator == 'random-parallel' )  THEN
2952       dist_nxl = nxl;  dist_nxr = nxr
2953       dist_nys = nys;  dist_nyn = nyn
2954       IF ( bc_lr == 'radiation/dirichlet' )  THEN
2955          dist_nxr    = MIN( nx - inflow_disturbance_begin, nxr )
2956          dist_nxl(1) = MAX( nx - inflow_disturbance_end, nxl )
2957       ELSEIF ( bc_lr == 'dirichlet/radiation' )  THEN
2958          dist_nxl    = MAX( inflow_disturbance_begin, nxl )
2959          dist_nxr(1) = MIN( inflow_disturbance_end, nxr )
2960       ELSEIF ( bc_lr == 'nested'  .OR.  bc_lr == 'nesting_offline' )  THEN
2961          dist_nxl    = MAX( inflow_disturbance_begin, nxl )
2962          dist_nxr    = MIN( nx - inflow_disturbance_begin, nxr )
2963       ENDIF
2964       IF ( bc_ns == 'dirichlet/radiation' )  THEN
2965          dist_nyn    = MIN( ny - inflow_disturbance_begin, nyn )
2966          dist_nys(1) = MAX( ny - inflow_disturbance_end, nys )
2967       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
2968          dist_nys    = MAX( inflow_disturbance_begin, nys )
2969          dist_nyn(1) = MIN( inflow_disturbance_end, nyn )
2970       ELSEIF ( bc_ns == 'nested'  .OR.  bc_ns == 'nesting_offline' )  THEN
2971          dist_nys    = MAX( inflow_disturbance_begin, nys )
2972          dist_nyn    = MIN( ny - inflow_disturbance_begin, nyn )
2973       ENDIF
2974    ELSE
2975       dist_nxl = 0;  dist_nxr = nx
2976       dist_nys = 0;  dist_nyn = ny
2977       IF ( bc_lr == 'radiation/dirichlet' )  THEN
2978          dist_nxr    = nx - inflow_disturbance_begin
2979          dist_nxl(1) = nx - inflow_disturbance_end
2980       ELSEIF ( bc_lr == 'dirichlet/radiation' )  THEN
2981          dist_nxl    = inflow_disturbance_begin
2982          dist_nxr(1) = inflow_disturbance_end
2983       ELSEIF ( bc_lr == 'nested'  .OR.  bc_lr == 'nesting_offline' )  THEN
2984          dist_nxr    = nx - inflow_disturbance_begin
2985          dist_nxl    = inflow_disturbance_begin
2986       ENDIF
2987       IF ( bc_ns == 'dirichlet/radiation' )  THEN
2988          dist_nyn    = ny - inflow_disturbance_begin
2989          dist_nys(1) = ny - inflow_disturbance_end
2990       ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
2991          dist_nys    = inflow_disturbance_begin
2992          dist_nyn(1) = inflow_disturbance_end
2993       ELSEIF ( bc_ns == 'nested'  .OR.  bc_ns == 'nesting_offline' )  THEN
2994          dist_nyn    = ny - inflow_disturbance_begin
2995          dist_nys    = inflow_disturbance_begin
2996       ENDIF
2997    ENDIF
2998
2999!
3000!-- A turbulent inflow requires Dirichlet conditions at the respective inflow
3001!-- boundary (so far, a turbulent inflow is realized from the left side only)
3002    IF ( turbulent_inflow  .AND.  bc_lr /= 'dirichlet/radiation' )  THEN
3003       message_string = 'turbulent_inflow = .T. requires a Dirichlet ' //      &
3004                        'condition at the inflow boundary'
3005       CALL message( 'check_parameters', 'PA0133', 1, 2, 0, 6, 0 )
3006    ENDIF
3007
3008!
3009!-- Turbulent inflow requires that 3d arrays have been cyclically filled with
3010!-- data from prerun in the first main run
3011    IF ( turbulent_inflow  .AND.  initializing_actions /= 'cyclic_fill'  .AND. &
3012         initializing_actions /= 'read_restart_data' )  THEN
3013       message_string = 'turbulent_inflow = .T. requires ' //                  &
3014                        'initializing_actions = ''cyclic_fill'' or ' //        &
3015                        'initializing_actions = ''read_restart_data'' '
3016       CALL message( 'check_parameters', 'PA0055', 1, 2, 0, 6, 0 )
3017    ENDIF
3018
3019!
3020!-- In case of turbulent inflow
3021    IF ( turbulent_inflow )  THEN
3022
3023!
3024!--    Calculate the index of the recycling plane
3025       IF ( recycling_width <= dx  .OR.  recycling_width >= nx * dx )  THEN
3026          WRITE( message_string, * )  'illegal value for recycling_width: ',   &
3027                                      recycling_width
3028          CALL message( 'check_parameters', 'PA0134', 1, 2, 0, 6, 0 )
3029       ENDIF
3030!
3031!--    Calculate the index
3032       recycling_plane = recycling_width / dx
3033!
3034!--   Check for correct input of recycling method for thermodynamic quantities
3035       IF ( TRIM( recycling_method_for_thermodynamic_quantities ) /= 'turbulent_fluctuation' .AND. &
3036            TRIM( recycling_method_for_thermodynamic_quantities ) /= 'absolute_value' )  THEN
3037          WRITE( message_string, * )  'unknown recycling method for thermodynamic quantities: ',   &
3038               TRIM( recycling_method_for_thermodynamic_quantities )
3039          CALL message( 'check_parameters', 'PA0184', 1, 2, 0, 6, 0 )
3040       ENDIF
3041
3042    ENDIF
3043
3044
3045    IF ( turbulent_outflow )  THEN
3046!
3047!--    Turbulent outflow requires Dirichlet conditions at the respective inflow
3048!--    boundary (so far, a turbulent outflow is realized at the right side only)
3049       IF ( bc_lr /= 'dirichlet/radiation' )  THEN
3050          message_string = 'turbulent_outflow = .T. requires ' //              &
3051                           'bc_lr = "dirichlet/radiation"'
3052          CALL message( 'check_parameters', 'PA0038', 1, 2, 0, 6, 0 )
3053       ENDIF
3054!
3055!--    The ouflow-source plane must lay inside the model domain
3056       IF ( outflow_source_plane < dx  .OR.  &
3057            outflow_source_plane > nx * dx )  THEN
3058          WRITE( message_string, * )  'illegal value for outflow_source'//     &
3059                                      '_plane: ', outflow_source_plane
3060          CALL message( 'check_parameters', 'PA0145', 1, 2, 0, 6, 0 )
3061       ENDIF
3062    ENDIF
3063
3064!
3065!-- Determine damping level index for 1D model
3066    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
3067       IF ( damp_level_1d == -1.0_wp )  THEN
3068          damp_level_1d     = zu(nzt+1)
3069          damp_level_ind_1d = nzt + 1
3070       ELSEIF ( damp_level_1d < 0.0_wp  .OR.  damp_level_1d > zu(nzt+1) )  THEN
3071          WRITE( message_string, * )  'damp_level_1d = ', damp_level_1d,       &
3072                 ' must be >= 0.0 and <= ', zu(nzt+1), '(zu(nzt+1))'
3073          CALL message( 'check_parameters', 'PA0136', 1, 2, 0, 6, 0 )
3074       ELSE
3075          DO  k = 1, nzt+1
3076             IF ( damp_level_1d <= zu(k) )  THEN
3077                damp_level_ind_1d = k
3078                EXIT
3079             ENDIF
3080          ENDDO
3081       ENDIF
3082    ENDIF
3083
3084!
3085!-- Check some other 1d-model parameters
3086    IF ( TRIM( mixing_length_1d ) /= 'as_in_3d_model'  .AND.                   &
3087         TRIM( mixing_length_1d ) /= 'blackadar' )  THEN
3088       message_string = 'mixing_length_1d = "' // TRIM( mixing_length_1d ) //  &
3089                        '" is unknown'
3090       CALL message( 'check_parameters', 'PA0137', 1, 2, 0, 6, 0 )
3091    ENDIF
3092    IF ( TRIM( dissipation_1d ) /= 'as_in_3d_model'  .AND.                     &
3093         TRIM( dissipation_1d ) /= 'detering'  .AND.                           &
3094         TRIM( dissipation_1d ) /= 'prognostic' )  THEN
3095       message_string = 'dissipation_1d = "' // TRIM( dissipation_1d ) //      &
3096                        '" is unknown'
3097       CALL message( 'check_parameters', 'PA0138', 1, 2, 0, 6, 0 )
3098    ENDIF
3099    IF ( TRIM( mixing_length_1d ) /= 'as_in_3d_model'  .AND.                   &
3100         TRIM( dissipation_1d ) == 'as_in_3d_model' )  THEN
3101       message_string = 'dissipation_1d = "' // TRIM( dissipation_1d ) //      &
3102                        '" requires mixing_length_1d = "as_in_3d_model"'
3103       CALL message( 'check_parameters', 'PA0485', 1, 2, 0, 6, 0 )
3104    ENDIF
3105
3106!
3107!-- Set time for the next user defined restart (time_restart is the
3108!-- internal parameter for steering restart events)
3109    IF ( restart_time /= 9999999.9_wp )  THEN
3110       IF ( restart_time > time_since_reference_point )  THEN
3111          time_restart = restart_time
3112       ENDIF
3113    ELSE
3114!
3115!--    In case of a restart run, set internal parameter to default (no restart)
3116!--    if the NAMELIST-parameter restart_time is omitted
3117       time_restart = 9999999.9_wp
3118    ENDIF
3119
3120!
3121!-- Check pressure gradient conditions
3122    IF ( dp_external  .AND.  conserve_volume_flow )  THEN
3123       WRITE( message_string, * )  'Both dp_external and conserve_volume_flo', &
3124            'w are .TRUE. but one of them must be .FALSE.'
3125       CALL message( 'check_parameters', 'PA0150', 1, 2, 0, 6, 0 )
3126    ENDIF
3127    IF ( dp_external )  THEN
3128       IF ( dp_level_b < zu(nzb)  .OR.  dp_level_b > zu(nzt) )  THEN
3129          WRITE( message_string, * )  'dp_level_b = ', dp_level_b, ' is out ', &
3130               ' of range [zu(nzb), zu(nzt)]'
3131          CALL message( 'check_parameters', 'PA0151', 1, 2, 0, 6, 0 )
3132       ENDIF
3133       IF ( .NOT. ANY( dpdxy /= 0.0_wp ) )  THEN
3134          WRITE( message_string, * )  'dp_external is .TRUE. but dpdxy is ze', &
3135               'ro, i.e. the external pressure gradient will not be applied'
3136          CALL message( 'check_parameters', 'PA0152', 0, 1, 0, 6, 0 )
3137       ENDIF
3138    ENDIF
3139    IF ( ANY( dpdxy /= 0.0_wp )  .AND.  .NOT.  dp_external )  THEN
3140       WRITE( message_string, * )  'dpdxy is nonzero but dp_external is ',     &
3141            '.FALSE., i.e. the external pressure gradient & will not be applied'
3142       CALL message( 'check_parameters', 'PA0153', 0, 1, 0, 6, 0 )
3143    ENDIF
3144    IF ( conserve_volume_flow )  THEN
3145       IF ( TRIM( conserve_volume_flow_mode ) == 'default' )  THEN
3146
3147          conserve_volume_flow_mode = 'initial_profiles'
3148
3149       ELSEIF ( TRIM( conserve_volume_flow_mode ) /= 'initial_profiles' .AND.  &
3150            TRIM( conserve_volume_flow_mode ) /= 'bulk_velocity' )  THEN
3151          WRITE( message_string, * )  'unknown conserve_volume_flow_mode: ',   &
3152               conserve_volume_flow_mode
3153          CALL message( 'check_parameters', 'PA0154', 1, 2, 0, 6, 0 )
3154       ENDIF
3155       IF ( (bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic')  .AND.                &
3156          TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
3157          WRITE( message_string, * )  'non-cyclic boundary conditions ',       &
3158               'require  conserve_volume_flow_mode = ''initial_profiles'''
3159          CALL message( 'check_parameters', 'PA0155', 1, 2, 0, 6, 0 )
3160       ENDIF
3161    ENDIF
3162    IF ( ( u_bulk /= 0.0_wp  .OR.  v_bulk /= 0.0_wp )  .AND.                   &
3163         ( .NOT. conserve_volume_flow  .OR.                                    &
3164         TRIM( conserve_volume_flow_mode ) /= 'bulk_velocity' ) )  THEN
3165       WRITE( message_string, * )  'nonzero bulk velocity requires ',          &
3166            'conserve_volume_flow = .T. and ',                                 &
3167            'conserve_volume_flow_mode = ''bulk_velocity'''
3168       CALL message( 'check_parameters', 'PA0157', 1, 2, 0, 6, 0 )
3169    ENDIF
3170
3171!
3172!-- Prevent empty time records in volume, cross-section and masked data in case
3173!-- of non-parallel netcdf-output in restart runs
3174    IF ( netcdf_data_format < 5 )  THEN
3175       IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
3176          do3d_time_count    = 0
3177          do2d_xy_time_count = 0
3178          do2d_xz_time_count = 0
3179          do2d_yz_time_count = 0
3180          domask_time_count  = 0
3181       ENDIF
3182    ENDIF
3183
3184
3185!
3186!-- Check roughness length, which has to be smaller than dz/2
3187    IF ( ( constant_flux_layer .OR.  &
3188           INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )       &
3189         .AND. roughness_length >= 0.5 * dz(1) )  THEN
3190       message_string = 'roughness_length must be smaller than dz/2'
3191       CALL message( 'check_parameters', 'PA0424', 1, 2, 0, 6, 0 )
3192    ENDIF
3193
3194#if defined( __parallel )
3195!
3196!-- Vertical nesting: check fine and coarse grid compatibility for data exchange
3197    IF ( vnested )  CALL vnest_check_parameters
3198#endif
3199
3200!
3201!-- Check if topography is read from file in case of complex terrain simulations
3202    IF ( complex_terrain  .AND.  TRIM( topography ) /= 'read_from_file' )  THEN
3203       message_string = 'complex_terrain requires topography' //               &
3204                        ' = ''read_from_file'''
3205       CALL message( 'check_parameters', 'PA0295', 1, 2, 0, 6, 0 )
3206    ENDIF
3207
3208!
3209!-- Check if vertical grid stretching is switched off in case of complex
3210!-- terrain simulations
3211    IF ( complex_terrain  .AND.                                                &
3212         ANY( dz_stretch_level_start /= -9999999.9_wp ) )  THEN
3213       message_string = 'Vertical grid stretching is not allowed for ' //      &
3214                        'complex_terrain = .T.'
3215       CALL message( 'check_parameters', 'PA0473', 1, 2, 0, 6, 0 )
3216    ENDIF
3217
3218    CALL location_message( 'checking parameters', 'finished' )
3219
3220 CONTAINS
3221
3222!------------------------------------------------------------------------------!
3223! Description:
3224! ------------
3225!> Check the length of data output intervals. In case of parallel NetCDF output
3226!> the time levels of the output files need to be fixed. Therefore setting the
3227!> output interval to 0.0s (usually used to output each timestep) is not
3228!> possible as long as a non-fixed timestep is used.
3229!------------------------------------------------------------------------------!
3230
3231    SUBROUTINE check_dt_do( dt_do, dt_do_name )
3232
3233       IMPLICIT NONE
3234
3235       CHARACTER (LEN=*), INTENT (IN) :: dt_do_name !< parin variable name
3236
3237       REAL(wp), INTENT (INOUT)       :: dt_do      !< data output interval
3238
3239       IF ( dt_do == 0.0_wp )  THEN
3240          IF ( dt_fixed )  THEN
3241             WRITE( message_string, '(A,F9.4,A)' )  'Output at every '  //     &
3242                    'timestep is wanted (' // dt_do_name // ' = 0.0).&'//      &
3243                    'The output interval is set to the fixed timestep dt '//   &
3244                    '= ', dt, 's.'
3245             CALL message( 'check_parameters', 'PA0060', 0, 0, 0, 6, 0 )
3246             dt_do = dt
3247          ELSE
3248             message_string = dt_do_name // ' = 0.0 while using a ' //         &
3249                              'variable timestep and parallel netCDF4 ' //     &
3250                              'is not allowed.'
3251             CALL message( 'check_parameters', 'PA0081', 1, 2, 0, 6, 0 )
3252          ENDIF
3253       ENDIF
3254
3255    END SUBROUTINE check_dt_do
3256
3257
3258
3259!------------------------------------------------------------------------------!
3260! Description:
3261! ------------
3262!> Set the bottom and top boundary conditions for humidity and scalars.
3263!------------------------------------------------------------------------------!
3264
3265    SUBROUTINE set_bc_scalars( sq, bc_b, bc_t, ibc_b, ibc_t, err_nr_b, err_nr_t )
3266
3267
3268       IMPLICIT NONE
3269
3270       CHARACTER (LEN=1)   ::  sq         !< name of scalar quantity
3271       CHARACTER (LEN=*)   ::  bc_b       !< bottom boundary condition
3272       CHARACTER (LEN=*)   ::  bc_t       !< top boundary condition
3273       CHARACTER (LEN=*)   ::  err_nr_b   !< error number if bottom bc is unknown
3274       CHARACTER (LEN=*)   ::  err_nr_t   !< error number if top bc is unknown
3275
3276       INTEGER(iwp)        ::  ibc_b      !< index for bottom boundary condition
3277       INTEGER(iwp)        ::  ibc_t      !< index for top boundary condition
3278
3279!
3280!--    Set Integer flags and check for possilbe errorneous settings for bottom
3281!--    boundary condition
3282       IF ( bc_b == 'dirichlet' )  THEN
3283          ibc_b = 0
3284       ELSEIF ( bc_b == 'neumann' )  THEN
3285          ibc_b = 1
3286       ELSE
3287          message_string = 'unknown boundary condition: bc_' // TRIM( sq ) //  &
3288                           '_b ="' // TRIM( bc_b ) // '"'
3289          CALL message( 'check_parameters', err_nr_b, 1, 2, 0, 6, 0 )
3290       ENDIF
3291!
3292!--    Set Integer flags and check for possilbe errorneous settings for top
3293!--    boundary condition
3294       IF ( bc_t == 'dirichlet' )  THEN
3295          ibc_t = 0
3296       ELSEIF ( bc_t == 'neumann' )  THEN
3297          ibc_t = 1
3298       ELSEIF ( bc_t == 'initial_gradient' )  THEN
3299          ibc_t = 2
3300       ELSEIF ( bc_t == 'nested'  .OR.  bc_t == 'nesting_offline' )  THEN
3301          ibc_t = 3
3302       ELSE
3303          message_string = 'unknown boundary condition: bc_' // TRIM( sq ) //  &
3304                           '_t ="' // TRIM( bc_t ) // '"'
3305          CALL message( 'check_parameters', err_nr_t, 1, 2, 0, 6, 0 )
3306       ENDIF
3307
3308
3309    END SUBROUTINE set_bc_scalars
3310
3311
3312
3313!------------------------------------------------------------------------------!
3314! Description:
3315! ------------
3316!> Check for consistent settings of bottom boundary conditions for humidity
3317!> and scalars.
3318!------------------------------------------------------------------------------!
3319
3320    SUBROUTINE check_bc_scalars( sq, bc_b, ibc_b,                      &
3321                                 err_nr_1, err_nr_2,                   &
3322                                 constant_flux, surface_initial_change )
3323
3324
3325       IMPLICIT NONE
3326
3327       CHARACTER (LEN=1)   ::  sq                       !< name of scalar quantity
3328       CHARACTER (LEN=*)   ::  bc_b                     !< bottom boundary condition
3329       CHARACTER (LEN=*)   ::  err_nr_1                 !< error number of first error
3330       CHARACTER (LEN=*)   ::  err_nr_2                 !< error number of second error
3331
3332       INTEGER(iwp)        ::  ibc_b                    !< index of bottom boundary condition
3333
3334       LOGICAL             ::  constant_flux            !< flag for constant-flux layer
3335
3336       REAL(wp)            ::  surface_initial_change   !< value of initial change at the surface
3337
3338!
3339!--    A given surface value implies Dirichlet boundary condition for
3340!--    the respective quantity. In this case specification of a constant flux is
3341!--    forbidden. However, an exception is made for large-scale forcing as well
3342!--    as land-surface model.
3343       IF ( .NOT. land_surface  .AND.  .NOT. large_scale_forcing )  THEN
3344          IF ( ibc_b == 0  .AND.  constant_flux )  THEN
3345             message_string = 'boundary condition: bc_' // TRIM( sq ) //       &
3346                              '_b ' // '= "' // TRIM( bc_b ) //                &
3347                              '" is not allowed with prescribed surface flux'
3348             CALL message( 'check_parameters', err_nr_1, 1, 2, 0, 6, 0 )
3349          ENDIF
3350       ENDIF
3351       IF ( constant_flux  .AND.  surface_initial_change /= 0.0_wp )  THEN
3352          WRITE( message_string, * )  'a prescribed surface flux is not allo', &
3353                 'wed with ', sq, '_surface_initial_change (/=0) = ',          &
3354                 surface_initial_change
3355          CALL message( 'check_parameters', err_nr_2, 1, 2, 0, 6, 0 )
3356       ENDIF
3357
3358
3359    END SUBROUTINE check_bc_scalars
3360
3361
3362
3363 END SUBROUTINE check_parameters
Note: See TracBrowser for help on using the repository browser.