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

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

Copyright updated to year 2021, interface pmc_sort removed to accelarate the nesting code

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