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

Last change on this file since 4511 was 4511, checked in by raasch, 4 years ago

chemistry decycling replaced by explicit setting of lateral boundary conditions

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