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

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

In case of turbulent inflow: Adjusted recycling_yshift, so that y-shift can be given as multiples of PE (as in y_shift for cyclic BCs).

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