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

Last change on this file since 4526 was 4514, checked in by suehring, 5 years ago

Bugfix in plant-canopy model for output of averaged transpiration rate after a restart; Revise check for output for plant heating rate and rename error message number; Surface-data output: enable output of mixing ratio and passive scalar concentration at the surface; Surface-data input: Add possibility to prescribe surface sensible and latent heat fluxes via static input file

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