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

Last change on this file since 4457 was 4444, checked in by raasch, 4 years ago

bugfix: cpp-directives for serial mode added

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