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

Last change on this file since 4336 was 4331, checked in by suehring, 5 years ago

New diagnostic output for 10-m wind speed; Diagnostic output of 2-m potential temperature moved to diagnostic output

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