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

Last change on this file since 4575 was 4565, checked in by oliver.maas, 4 years ago

added new surface temperature forcing method for bc_pt_b = 'dirichlet': surface temperature pt(0) can be linearly increased by pt_surface_heating_rate (in K/h)

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