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

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

Deleted parameter recycling_yshift. y-shift in case of non-cyclic boundary conditions and turbulent_inflow = .TRUE. is now steered by parameter y_shift, that is also used in case of cyclic boundary conditions.

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