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

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