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

Last change on this file since 112 was 110, checked in by raasch, 17 years ago

New:
---
Allows runs for a coupled atmosphere-ocean LES,
coupling frequency is controlled by new d3par-parameter dt_coupling,
the coupling mode (atmosphere_to_ocean or ocean_to_atmosphere) for the
respective processes is read from environment variable coupling_mode,
which is set by the mpiexec-command,
communication between the two models is done using the intercommunicator
comm_inter,
local files opened by the ocean model get the additional suffic "_O".
Assume saturation at k=nzb_s_inner(j,i) for atmosphere coupled to ocean.

A momentum flux can be set as top boundary condition using the new
inipar parameter top_momentumflux_u|v.

Non-cyclic boundary conditions can be used along all horizontal directions.

Quantities w*p* and w"e can be output as vertical profiles.

Initial profiles are reset to constant profiles in case that initializing_actions /= 'set_constant_profiles'. (init_rankine)

Optionally calculate km and kh from initial TKE e_init.

Changed:


Remaining variables iran changed to iran_part (advec_particles, init_particles).

In case that the presure solver is not called for every Runge-Kutta substep
(call_psolver_at_all_substeps = .F.), it is called after the first substep
instead of the last. In that case, random perturbations are also added to the
velocity field after the first substep.

Initialization of km,kh = 0.00001 for ocean = .T. (for ocean = .F. it remains 0.01).

Allow data_output_pr= q, wq, w"q", w*q* for humidity = .T. (instead of cloud_physics = .T.).

Errors:


Bugs from code parts for non-cyclic boundary conditions are removed: loops for
u and v are starting from index nxlu, nysv, respectively. The radiation boundary
condition is used for every Runge-Kutta substep. Velocity phase speeds for
the radiation boundary conditions are calculated for the first Runge-Kutta
substep only and reused for the further substeps. New arrays c_u, c_v, and c_w
are defined for this purpose. Several index errors are removed from the
radiation boundary condition code parts. Upper bounds for calculating
u_0 and v_0 (in production_e) are nxr+1 and nyn+1 because otherwise these
values are not available in case of non-cyclic boundary conditions.

+dots_num_palm in module user, +module netcdf_control in user_init (both in user_interface)

Bugfix: wrong sign removed from the buoyancy production term in the case use_reference = .T. (production_e)

Bugfix: Error message concerning output of particle concentration (pc) modified (check_parameters).

Bugfix: Rayleigh damping for ocean fixed.

  • Property svn:keywords set to Id
File size: 109.0 KB
Line 
1 SUBROUTINE check_parameters
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6!
7!
8! Former revisions:
9! -----------------
10! $Id: check_parameters.f90 110 2007-10-05 05:13:14Z raasch $
11!
12! 109 2007-08-28 15:26:47Z letzel
13! Check coupling_mode and set default (obligatory) values (like boundary
14! conditions for temperature and fluxes) in case of coupled runs.
15! +profiles for w*p* and w"e
16! Bugfix: Error message concerning output of particle concentration (pc)
17! modified
18! More checks and more default values for coupled runs
19! allow data_output_pr= q, wq, w"q", w*q* for humidity = .T. (instead of
20! cloud_physics = .T.)
21! Rayleigh damping for ocean fixed.
22! Check and, if necessary, set default value for dt_coupling
23!
24! 97 2007-06-21 08:23:15Z raasch
25! Initial salinity profile is calculated, salinity boundary conditions are
26! checked,
27! z_max_do1d is checked only in case of ocean = .f.,
28! +initial temperature and geostrophic velocity profiles for the ocean version,
29! use_pt_reference renamed use_reference
30!
31! 89 2007-05-25 12:08:31Z raasch
32! Check for user-defined profiles
33!
34! 75 2007-03-22 09:54:05Z raasch
35! "by_user" allowed as initializing action, -data_output_ts,
36! leapfrog with non-flat topography not allowed any more, loop_optimization
37! and pt_reference are checked, moisture renamed humidity,
38! output of precipitation amount/rate and roughnes length + check
39! possible negative humidities are avoided in initial profile,
40! dirichlet/neumann changed to dirichlet/radiation, etc.,
41! revision added to run_description_header
42!
43! 20 2007-02-26 00:12:32Z raasch
44! Temperature and humidity gradients at top are now calculated for nzt+1,
45! top_heatflux and respective boundary condition bc_pt_t is checked
46!
47! RCS Log replace by Id keyword, revision history cleaned up
48!
49! Revision 1.61  2006/08/04 14:20:25  raasch
50! do2d_unit and do3d_unit now defined as 2d-arrays, check of
51! use_upstream_for_tke, default value for dt_dopts,
52! generation of file header moved from routines palm and header to here
53!
54! Revision 1.1  1997/08/26 06:29:23  raasch
55! Initial revision
56!
57!
58! Description:
59! ------------
60! Check control parameters and deduce further quantities.
61!------------------------------------------------------------------------------!
62
63    USE arrays_3d
64    USE constants
65    USE control_parameters
66    USE grid_variables
67    USE indices
68    USE model_1d
69    USE netcdf_control
70    USE particle_attributes
71    USE pegrid
72    USE profil_parameter
73    USE statistics
74    USE transpose_indices
75
76    IMPLICIT NONE
77
78    CHARACTER (LEN=1)   ::  sq
79    CHARACTER (LEN=6)   ::  var
80    CHARACTER (LEN=7)   ::  unit
81    CHARACTER (LEN=8)   ::  date
82    CHARACTER (LEN=10)  ::  time
83    CHARACTER (LEN=40)  ::  coupling_string
84    CHARACTER (LEN=100) ::  action
85
86    INTEGER ::  i, ilen, intervals, iremote = 0, iter, j, k, nnxh, nnyh, &
87         position, prec
88    LOGICAL ::  found, ldum
89    REAL    ::  gradient, maxn, maxp, remote = 0.0
90
91!
92!-- Warning, if host is not set
93    IF ( host(1:1) == ' ' )  THEN
94       IF ( myid == 0 )  THEN
95          PRINT*, '+++ WARNING: check_parameters:'
96          PRINT*, '    "host" is not set. Please check that environment', &
97                       ' variable "localhost"'
98          PRINT*, '    is set before running PALM'
99       ENDIF
100    ENDIF
101
102!
103!-- Check the coupling mode
104    IF ( coupling_mode /= 'uncoupled'            .AND.  &
105         coupling_mode /= 'atmosphere_to_ocean'  .AND.  &
106         coupling_mode /= 'ocean_to_atmosphere' )  THEN
107       IF ( myid == 0 )  THEN
108          PRINT*, '+++ check_parameters:'
109          PRINT*, '    illegal coupling mode: ', TRIM( coupling_mode )
110       ENDIF
111       CALL local_stop
112    ENDIF
113
114!
115!-- Check dt_coupling, restart_time, dt_restart, end_time, dx, dy, nx and ny
116    IF ( coupling_mode /= 'uncoupled' )  THEN
117       IF ( dt_coupling == 9999999.9 )  THEN
118          IF ( myid == 0 )  THEN
119             PRINT*, '+++ check_parameters:'
120             PRINT*, '    dt_coupling is not set but required for coupling ', &
121                  'mode: ', TRIM( coupling_mode )
122          ENDIF
123          CALL local_stop
124       ENDIF
125#if defined( __parallel )  &&  defined( __mpi2 )
126       CALL MPI_SEND( dt_coupling, 1, MPI_REAL, myid, 11, comm_inter, ierr )
127       CALL MPI_RECV( remote, 1, MPI_REAL, myid, 11, comm_inter, status, ierr )
128       IF ( dt_coupling /= remote )  THEN
129          IF ( myid == 0 )  THEN
130             PRINT*, '+++ check_parameters:'
131             PRINT*, '    TRIM( coupling_mode ): dt_coupling = ', dt_coupling
132             PRINT*, '    is not equal to dt_coupling_remote = ', remote
133          ENDIF
134          CALL local_stop
135       ENDIF
136       IF ( dt_coupling <= 0.0 )  THEN
137          CALL MPI_SEND( dt_max, 1, MPI_REAL, myid, 19, comm_inter, ierr )
138          CALL MPI_RECV( remote, 1, MPI_REAL, myid, 19, comm_inter, status, &
139               ierr )
140          dt_coupling = MAX( dt_max, remote )
141          IF ( myid == 0 )  THEN
142             PRINT*, '+++ check_parameters:'
143             PRINT*, '    TRIM( coupling_mode ): dt_coupling <= 0.0'
144             PRINT*, '    is not allowed and is reset to MAX(dt_max(A,O)) = ', &
145                  dt_coupling
146          ENDIF
147       ENDIF
148       CALL MPI_SEND( restart_time, 1, MPI_REAL, myid, 12, comm_inter, ierr )
149       CALL MPI_RECV( remote, 1, MPI_REAL, myid, 12, comm_inter, status, ierr )
150       IF ( restart_time /= remote )  THEN
151          IF ( myid == 0 )  THEN
152             PRINT*, '+++ check_parameters:'
153             PRINT*, '    TRIM( coupling_mode ): restart_time = ', restart_time
154             PRINT*, '    is not equal to restart_time_remote = ', remote
155          ENDIF
156          CALL local_stop
157       ENDIF
158       CALL MPI_SEND( dt_restart, 1, MPI_REAL, myid, 13, comm_inter, ierr )
159       CALL MPI_RECV( remote, 1, MPI_REAL, myid, 13, comm_inter, status, ierr )
160       IF ( dt_restart /= remote )  THEN
161          IF ( myid == 0 )  THEN
162             PRINT*, '+++ check_parameters:'
163             PRINT*, '    TRIM( coupling_mode ): dt_restart = ', dt_restart
164             PRINT*, '    is not equal to dt_restart_remote = ', remote
165          ENDIF
166          CALL local_stop
167       ENDIF
168       CALL MPI_SEND( end_time, 1, MPI_REAL, myid, 14, comm_inter, ierr )
169       CALL MPI_RECV( remote, 1, MPI_REAL, myid, 14, comm_inter, status, ierr )
170       IF ( end_time /= remote )  THEN
171          IF ( myid == 0 )  THEN
172             PRINT*, '+++ check_parameters:'
173             PRINT*, '    TRIM( coupling_mode ): end_time = ', end_time
174             PRINT*, '    is not equal to end_time_remote = ', remote
175          ENDIF
176          CALL local_stop
177       ENDIF
178       CALL MPI_SEND( dx, 1, MPI_REAL, myid, 15, comm_inter, ierr )
179       CALL MPI_RECV( remote, 1, MPI_REAL, myid, 15, comm_inter, status, ierr )
180       IF ( dx /= remote )  THEN
181          IF ( myid == 0 )  THEN
182             PRINT*, '+++ check_parameters:'
183             PRINT*, '    TRIM( coupling_mode ): dx = ', dx
184             PRINT*, '    is not equal to dx_remote = ', remote
185          ENDIF
186          CALL local_stop
187       ENDIF
188       CALL MPI_SEND( dy, 1, MPI_REAL, myid, 16, comm_inter, ierr )
189       CALL MPI_RECV( remote, 1, MPI_REAL, myid, 16, comm_inter, status, ierr )
190       IF ( dy /= remote )  THEN
191          IF ( myid == 0 )  THEN
192             PRINT*, '+++ check_parameters:'
193             PRINT*, '    TRIM( coupling_mode ): dy = ', dy
194             PRINT*, '    is not equal to dy_remote = ', remote
195          ENDIF
196          CALL local_stop
197       ENDIF
198       CALL MPI_SEND( nx, 1, MPI_INTEGER, myid, 17, comm_inter, ierr )
199       CALL MPI_RECV( iremote, 1, MPI_INTEGER, myid, 17, comm_inter, status, &
200            ierr )
201       IF ( nx /= iremote )  THEN
202          IF ( myid == 0 )  THEN
203             PRINT*, '+++ check_parameters:'
204             PRINT*, '    TRIM( coupling_mode ): nx = ', nx
205             PRINT*, '    is not equal to nx_remote = ', iremote
206          ENDIF
207          CALL local_stop
208       ENDIF
209       CALL MPI_SEND( ny, 1, MPI_INTEGER, myid, 18, comm_inter, ierr )
210       CALL MPI_RECV( iremote, 1, MPI_INTEGER, myid, 18, comm_inter, status, &
211            ierr )
212       IF ( ny /= iremote )  THEN
213          IF ( myid == 0 )  THEN
214             PRINT*, '+++ check_parameters:'
215             PRINT*, '    TRIM( coupling_mode ): ny = ', ny
216             PRINT*, '    is not equal to ny_remote = ', iremote
217          ENDIF
218          CALL local_stop
219       ENDIF
220#endif
221    ENDIF
222
223#if defined( __parallel )  &&  defined( __mpi2 )
224!
225!-- Exchange via intercommunicator
226    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
227       CALL MPI_SEND( humidity, &
228            1, MPI_LOGICAL, myid, 19, comm_inter, ierr )
229    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
230       CALL MPI_RECV( humidity_remote, &
231            1, MPI_LOGICAL, myid, 19, comm_inter, status, ierr )
232    ENDIF
233#endif
234
235
236!
237!-- Generate the file header which is used as a header for most of PALM's
238!-- output files
239    CALL DATE_AND_TIME( date, time )
240    run_date = date(7:8)//'-'//date(5:6)//'-'//date(3:4)
241    run_time = time(1:2)//':'//time(3:4)//':'//time(5:6)
242    IF ( coupling_mode == 'uncoupled' )  THEN
243       coupling_string = ''
244    ELSEIF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
245       coupling_string = ' coupled (atmosphere)'
246    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
247       coupling_string = ' coupled (ocean)'
248    ENDIF       
249
250    WRITE ( run_description_header,                                        &
251                             '(A,2X,A,2X,A,A,A,I2.2,A,2X,A,A,2X,A,1X,A)' ) &
252              TRIM( version ), TRIM( revision ), 'run: ',                  &
253              TRIM( run_identifier ), '.', runnr, TRIM( coupling_string ), &
254              'host: ', TRIM( host ), run_date, run_time
255
256!
257!-- Check the general loop optimization method
258    IF ( loop_optimization == 'default' )  THEN
259       IF ( host(1:3) == 'nec' )  THEN
260          loop_optimization = 'vector'
261       ELSE
262          loop_optimization = 'cache'
263       ENDIF
264    ENDIF
265    IF ( loop_optimization /= 'noopt'  .AND.  loop_optimization /= 'cache' &
266         .AND.  loop_optimization /= 'vector' )  THEN
267       IF ( myid == 0 )  THEN
268          PRINT*, '+++ check_parameters:'
269          PRINT*, '    illegal value given for loop_optimization: ', &
270                  TRIM( loop_optimization )
271       ENDIF
272       CALL local_stop
273    ENDIF
274
275!
276!-- Check topography setting (check for illegal parameter combinations)
277    IF ( topography /= 'flat' )  THEN
278       action = ' '
279       IF ( scalar_advec /= 'pw-scheme' )  THEN
280          WRITE( action, '(A,A)' )  'scalar_advec = ', scalar_advec
281       ENDIF
282       IF ( momentum_advec /= 'pw-scheme' )  THEN
283          WRITE( action, '(A,A)' )  'momentum_advec = ', momentum_advec
284       ENDIF
285       IF ( timestep_scheme(1:8) == 'leapfrog' )  THEN
286          WRITE( action, '(A,A)' )  'timestep_scheme = ', timestep_scheme
287       ENDIF
288       IF ( psolver == 'multigrid'  .OR.  psolver == 'sor' )  THEN
289          WRITE( action, '(A,A)' )  'psolver = ', psolver
290       ENDIF
291       IF ( sloping_surface )  THEN
292          WRITE( action, '(A)' )  'sloping surface = .TRUE.'
293       ENDIF
294       IF ( galilei_transformation )  THEN
295          WRITE( action, '(A)' )  'galilei_transformation = .TRUE.'
296       ENDIF
297       IF ( cloud_physics )  THEN
298          WRITE( action, '(A)' )  'cloud_physics = .TRUE.'
299       ENDIF
300       IF ( cloud_droplets )  THEN
301          WRITE( action, '(A)' )  'cloud_droplets = .TRUE.'
302       ENDIF
303       IF ( humidity )  THEN
304          WRITE( action, '(A)' )  'humidity = .TRUE.'
305       ENDIF
306       IF ( .NOT. prandtl_layer )  THEN
307          WRITE( action, '(A)' )  'prandtl_layer = .FALSE.'
308       ENDIF
309       IF ( action /= ' ' )  THEN
310          IF ( myid == 0 )  THEN
311             PRINT*, '+++ check_parameters:'
312             PRINT*, '    a non-flat topography does not allow ', TRIM( action )
313          ENDIF
314          CALL local_stop
315       ENDIF
316    ENDIF
317
318!
319!-- Check ocean setting
320    IF ( ocean )  THEN
321       action = ' '
322       IF ( timestep_scheme(1:8) == 'leapfrog' )  THEN
323          WRITE( action, '(A,A)' )  'timestep_scheme = ', timestep_scheme
324       ENDIF
325       IF ( momentum_advec == 'ups-scheme' )  THEN
326          WRITE( action, '(A,A)' )  'momentum_advec = ', momentum_advec
327       ENDIF
328       IF ( action /= ' ' )  THEN
329          IF ( myid == 0 )  THEN
330             PRINT*, '+++ check_parameters:'
331             PRINT*, '    ocean = .T. does not allow ', TRIM( action )
332          ENDIF
333          CALL local_stop
334       ENDIF
335    ENDIF
336
337!
338!-- Check whether there are any illegal values
339!-- Pressure solver:
340    IF ( psolver /= 'poisfft'  .AND.  psolver /= 'poisfft_hybrid'  .AND. &
341         psolver /= 'sor'  .AND.  psolver /= 'multigrid' )  THEN
342       IF ( myid == 0 )  THEN
343          PRINT*, '+++ check_parameters:'
344          PRINT*, '    unknown solver for perturbation pressure: psolver=', &
345                  psolver
346       ENDIF
347       CALL local_stop
348    ENDIF
349
350#if defined( __parallel )
351    IF ( psolver == 'poisfft_hybrid'  .AND.  pdims(2) /= 1 )  THEN
352       IF ( myid == 0 )  THEN
353          PRINT*, '+++ check_parameters:'
354          PRINT*, '    psolver="', TRIM( psolver ), '" only works for a ', &
355                       '1d domain-decomposition along x'
356          PRINT*, '    please do not set npey/=1 in the parameter file'
357       ENDIF
358       CALL local_stop
359    ENDIF
360    IF ( ( psolver == 'poisfft_hybrid'  .OR.  psolver == 'multigrid' )  .AND.  &
361         ( nxra > nxr  .OR.  nyna > nyn  .OR.  nza > nz ) )  THEN
362       IF ( myid == 0 )  THEN
363          PRINT*, '+++ check_parameters:'
364          PRINT*, '    psolver="', TRIM( psolver ), '" does not work for ', &
365                       'subdomains with unequal size'
366          PRINT*, '    please set grid_matching = ''strict'' in the parameter',&
367                       ' file'
368       ENDIF
369       CALL local_stop
370    ENDIF
371#else
372    IF ( psolver == 'poisfft_hybrid' )  THEN
373       IF ( myid == 0 )  THEN
374          PRINT*, '+++ check_parameters:'
375          PRINT*, '    psolver="', TRIM( psolver ), '" only works for a ', &
376                       'parallel environment'
377       ENDIF
378       CALL local_stop
379    ENDIF
380#endif
381
382    IF ( psolver == 'multigrid' )  THEN
383       IF ( cycle_mg == 'w' )  THEN
384          gamma_mg = 2
385       ELSEIF ( cycle_mg == 'v' )  THEN
386          gamma_mg = 1
387       ELSE
388          IF ( myid == 0 )  THEN
389             PRINT*, '+++ check_parameters:'
390             PRINT*, '    unknown multigrid cycle: cycle_mg=', cycle_mg
391          ENDIF
392          CALL local_stop
393       ENDIF
394    ENDIF
395
396    IF ( fft_method /= 'singleton-algorithm'  .AND.  &
397         fft_method /= 'temperton-algorithm'  .AND.  &
398         fft_method /= 'system-specific' )  THEN
399       IF ( myid == 0 )  THEN
400          PRINT*, '+++ check_parameters:'
401          PRINT*, '    unknown fft-algorithm: fft_method=', fft_method
402       ENDIF
403       CALL local_stop
404    ENDIF
405
406!
407!-- Advection schemes:
408    IF ( momentum_advec /= 'pw-scheme' .AND. momentum_advec /= 'ups-scheme' ) &
409    THEN
410       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  unknown advection ', &
411                                 'scheme: momentum_advec=', momentum_advec
412       CALL local_stop
413    ENDIF
414    IF ( ( momentum_advec == 'ups-scheme'  .OR.  scalar_advec == 'ups-scheme' )&
415                                      .AND.  timestep_scheme /= 'euler' )  THEN
416       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  momentum_advec=', &
417                                 momentum_advec, ' is not allowed with ', &
418                                 'timestep_scheme=', timestep_scheme
419       CALL local_stop
420    ENDIF
421
422    IF ( scalar_advec /= 'pw-scheme'  .AND.  scalar_advec /= 'bc-scheme'  .AND.&
423         scalar_advec /= 'ups-scheme' )  THEN
424       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  unknown advection ', &
425                                 'scheme: scalar_advec=', scalar_advec
426       CALL local_stop
427    ENDIF
428
429    IF ( use_sgs_for_particles  .AND.  .NOT. use_upstream_for_tke )  THEN
430       use_upstream_for_tke = .TRUE.
431       IF ( myid == 0 )  THEN
432          PRINT*, '+++ WARNING: check_parameters:  use_upstream_for_tke set ', &
433                       '.TRUE. because use_sgs_for_particles = .TRUE.'
434       ENDIF
435    ENDIF
436
437    IF ( use_upstream_for_tke  .AND.  timestep_scheme(1:8) == 'leapfrog' )  THEN
438       IF ( myid == 0 )  THEN
439          PRINT*, '+++ check_parameters:  use_upstream_for_tke = .TRUE. ', &
440                       'not allowed with timestep_scheme=', timestep_scheme
441       ENDIF
442       CALL local_stop
443    ENDIF
444
445!
446!-- Timestep schemes:
447    SELECT CASE ( TRIM( timestep_scheme ) )
448
449       CASE ( 'euler' )
450          intermediate_timestep_count_max = 1
451          asselin_filter_factor           = 0.0
452
453       CASE ( 'leapfrog', 'leapfrog+euler' )
454          intermediate_timestep_count_max = 1
455
456       CASE ( 'runge-kutta-2' )
457          intermediate_timestep_count_max = 2
458          asselin_filter_factor           = 0.0
459
460       CASE ( 'runge-kutta-3' )
461          intermediate_timestep_count_max = 3
462          asselin_filter_factor           = 0.0
463
464       CASE DEFAULT
465          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  unknown timestep ',&
466                                    'scheme: timestep_scheme=', timestep_scheme
467          CALL local_stop
468
469    END SELECT
470
471    IF ( scalar_advec == 'ups-scheme'  .AND.  timestep_scheme(1:5) == 'runge' )&
472    THEN
473       IF ( myid == 0 )  THEN
474          PRINT*, '+++ check_parameters:  scalar advection scheme "', &
475                                          TRIM( scalar_advec ), '"'
476          PRINT*, '    does not work with timestep_scheme "', &
477                                          TRIM( timestep_scheme ), '"'
478       ENDIF
479       CALL local_stop
480    ENDIF
481
482    IF ( momentum_advec /= 'pw-scheme' .AND. timestep_scheme(1:5) == 'runge' ) &
483    THEN
484       IF ( myid == 0 )  THEN
485          PRINT*, '+++ check_parameters:  momentum advection scheme "', &
486                                          TRIM( momentum_advec ), '"'
487          PRINT*, '    does not work with timestep_scheme "', &
488                                          TRIM( timestep_scheme ), '"'
489       ENDIF
490       CALL local_stop
491    ENDIF
492
493
494    IF ( initializing_actions == ' ' )  THEN
495       IF ( myid == 0 )  THEN
496          PRINT*, '+++ check parameters:'
497          PRINT*, '    no value found for initializing_actions'
498       ENDIF
499       CALL local_stop
500    ENDIF
501
502    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
503!
504!--    No model continuation run; several initialising actions are possible
505       action = initializing_actions
506       DO WHILE ( TRIM( action ) /= '' )
507          position = INDEX( action, ' ' )
508          SELECT CASE ( action(1:position-1) )
509
510             CASE ( 'set_constant_profiles', 'set_1d-model_profiles', &
511                    'by_user', 'initialize_vortex',     'initialize_ptanom' )
512                action = action(position+1:)
513
514             CASE DEFAULT
515                IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializi', &
516                               'ng_action unkown or not allowed: action = "', &
517                               TRIM(action), '"'
518                CALL local_stop
519
520          END SELECT
521       ENDDO
522    ENDIF
523    IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0  .AND. &
524         INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
525       IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializing_actions', &
526          '"set_constant_profiles" and "set_1d-model_profiles" are not', &
527          ' allowed simultaneously'
528       CALL local_stop
529    ENDIF
530    IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0  .AND. &
531         INDEX( initializing_actions, 'by_user' ) /= 0 )  THEN
532       IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializing_actions', &
533          '"set_constant_profiles" and "by_user" are not', &
534          ' allowed simultaneously'
535       CALL local_stop
536    ENDIF
537    IF ( INDEX( initializing_actions, 'by_user' ) /= 0  .AND. &
538         INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
539       IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializing_actions', &
540          '"by_user" and "set_1d-model_profiles" are not', &
541          ' allowed simultaneously'
542       CALL local_stop
543    ENDIF
544
545    IF ( cloud_physics  .AND.  .NOT. humidity )  THEN
546       IF ( myid == 0 )  PRINT*, '+++ check_parameters: cloud_physics =', &
547                                 cloud_physics, ' is not allowed with ',  &
548                                 'humidity =', humidity
549       CALL local_stop
550    ENDIF
551
552    IF ( precipitation  .AND.  .NOT.  cloud_physics )  THEN
553       IF ( myid == 0 )  PRINT*, '+++ check_parameters: precipitation =', &
554                                 precipitation, ' is not allowed with ',  &
555                                 'cloud_physics =', cloud_physics
556       CALL local_stop
557    ENDIF
558
559    IF ( humidity  .AND.  sloping_surface )  THEN
560       IF ( myid == 0 )  PRINT*, '+++ check_parameters: humidity = TRUE', &
561                                 'and hang = TRUE are not',               &
562                                 ' allowed simultaneously' 
563       CALL local_stop       
564    ENDIF
565
566    IF ( humidity  .AND.  scalar_advec == 'ups-scheme' )  THEN
567       IF ( myid == 0 )  PRINT*, '+++ check_parameters: UPS-scheme', &
568                                 'is not implemented for humidity'
569       CALL local_stop       
570    ENDIF
571
572    IF ( passive_scalar  .AND.  humidity )  THEN
573       IF ( myid == 0 )  PRINT*, '+++ check_parameters: humidity = TRUE and', &
574                                 'passive_scalar = TRUE is not allowed ',     &
575                                 'simultaneously'
576       CALL local_stop
577    ENDIF
578
579    IF ( passive_scalar  .AND.  scalar_advec == 'ups-scheme' )  THEN
580       IF ( myid == 0 )  PRINT*, '+++ check_parameters: UPS-scheme', &
581                                 'is not implemented for passive_scalar'
582       CALL local_stop       
583    ENDIF
584
585    IF ( grid_matching /= 'strict'  .AND.  grid_matching /= 'match' )  THEN
586       IF ( myid == 0 )  PRINT*, '+++ check_parameters: illegal value "', &
587                                 TRIM( grid_matching ),                   &
588                                 '" found for parameter grid_matching'
589       CALL local_stop       
590    ENDIF
591
592!
593!-- In case of no model continuation run, check initialising parameters and
594!-- deduce further quantities
595    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
596
597!
598!--    Initial profiles for 1D and 3D model, respectively
599       u_init  = ug_surface
600       v_init  = vg_surface
601       pt_init = pt_surface
602       IF ( humidity )        q_init  = q_surface
603       IF ( ocean )           sa_init = sa_surface
604       IF ( passive_scalar )  q_init  = s_surface
605
606!
607!--
608!--    If required, compute initial profile of the geostrophic wind
609!--    (component ug)
610       i = 1
611       gradient = 0.0
612
613       IF ( .NOT. ocean )  THEN
614
615          ug_vertical_gradient_level_ind(1) = 0
616          ug(0) = ug_surface
617          DO  k = 1, nzt+1
618             IF ( ug_vertical_gradient_level(i) < zu(k)  .AND. &
619                  ug_vertical_gradient_level(i) >= 0.0 )  THEN
620                gradient = ug_vertical_gradient(i) / 100.0
621                ug_vertical_gradient_level_ind(i) = k - 1
622                i = i + 1
623                IF ( i > 10 )  THEN
624                   IF ( myid == 0 )  THEN
625                      PRINT*, '+++ check_parameters: upper bound 10 of array', &
626                              ' "ug_vertical_gradient_level_ind" exceeded'
627                   ENDIF
628                   CALL local_stop
629                ENDIF
630             ENDIF
631             IF ( gradient /= 0.0 )  THEN
632                IF ( k /= 1 )  THEN
633                   ug(k) = ug(k-1) + dzu(k) * gradient
634                ELSE
635                   ug(k) = ug_surface + 0.5 * dzu(k) * gradient
636                ENDIF
637             ELSE
638                ug(k) = ug(k-1)
639             ENDIF
640          ENDDO
641
642       ELSE
643
644          ug_vertical_gradient_level_ind(1) = nzt+1
645          DO  k = nzt, 0, -1
646             IF ( ug_vertical_gradient_level(i) > zu(k)  .AND. &
647                  ug_vertical_gradient_level(i) <= 0.0 )  THEN
648                gradient = ug_vertical_gradient(i) / 100.0
649                ug_vertical_gradient_level_ind(i) = k + 1
650                i = i + 1
651                IF ( i > 10 )  THEN
652                   IF ( myid == 0 )  THEN
653                      PRINT*, '+++ check_parameters: upper bound 10 of array', &
654                              ' "ug_vertical_gradient_level_ind" exceeded'
655                   ENDIF
656                   CALL local_stop
657                ENDIF
658             ENDIF
659             IF ( gradient /= 0.0 )  THEN
660                IF ( k /= nzt )  THEN
661                   ug(k) = ug(k+1) - dzu(k+1) * gradient
662                ELSE
663                   ug(k)   = ug_surface - 0.5 * dzu(k+1) * gradient
664                   ug(k+1) = ug_surface + 0.5 * dzu(k+1) * gradient
665                ENDIF
666             ELSE
667                ug(k) = ug(k+1)
668             ENDIF
669          ENDDO
670
671       ENDIF
672
673       u_init = ug
674
675!
676!--    In case of no given gradients for ug, choose a vanishing gradient
677       IF ( ug_vertical_gradient_level(1) == -9999999.9 )  THEN
678          ug_vertical_gradient_level(1) = 0.0
679       ENDIF 
680
681!
682!--
683!--    If required, compute initial profile of the geostrophic wind
684!--    (component vg)
685       i = 1
686       gradient = 0.0
687
688       IF ( .NOT. ocean )  THEN
689
690          vg_vertical_gradient_level_ind(1) = 0
691          vg(0) = vg_surface
692          DO  k = 1, nzt+1
693             IF ( vg_vertical_gradient_level(i) < zu(k)  .AND. &
694                  vg_vertical_gradient_level(i) >= 0.0 )  THEN
695                gradient = vg_vertical_gradient(i) / 100.0
696                vg_vertical_gradient_level_ind(i) = k - 1
697                i = i + 1
698                IF ( i > 10 )  THEN
699                   IF ( myid == 0 )  THEN
700                      PRINT*, '+++ check_parameters: upper bound 10 of array', &
701                              ' "vg_vertical_gradient_level_ind" exceeded'
702                   ENDIF
703                   CALL local_stop
704                ENDIF
705             ENDIF
706             IF ( gradient /= 0.0 )  THEN
707                IF ( k /= 1 )  THEN
708                   vg(k) = vg(k-1) + dzu(k) * gradient
709                ELSE
710                   vg(k) = vg_surface + 0.5 * dzu(k) * gradient
711                ENDIF
712             ELSE
713                vg(k) = vg(k-1)
714             ENDIF
715          ENDDO
716
717       ELSE
718
719          vg_vertical_gradient_level_ind(1) = 0
720          DO  k = nzt, 0, -1
721             IF ( vg_vertical_gradient_level(i) > zu(k)  .AND. &
722                  vg_vertical_gradient_level(i) <= 0.0 )  THEN
723                gradient = vg_vertical_gradient(i) / 100.0
724                vg_vertical_gradient_level_ind(i) = k + 1
725                i = i + 1
726                IF ( i > 10 )  THEN
727                   IF ( myid == 0 )  THEN
728                      PRINT*, '+++ check_parameters: upper bound 10 of array', &
729                              ' "vg_vertical_gradient_level_ind" exceeded'
730                   ENDIF
731                   CALL local_stop
732                ENDIF
733             ENDIF
734             IF ( gradient /= 0.0 )  THEN
735                IF ( k /= nzt )  THEN
736                   vg(k) = vg(k+1) - dzu(k+1) * gradient
737                ELSE
738                   vg(k)   = vg_surface - 0.5 * dzu(k+1) * gradient
739                   vg(k+1) = vg_surface + 0.5 * dzu(k+1) * gradient
740                ENDIF
741             ELSE
742                vg(k) = vg(k+1)
743             ENDIF
744          ENDDO
745
746       ENDIF
747
748       v_init = vg
749 
750!
751!--    In case of no given gradients for vg, choose a vanishing gradient
752       IF ( vg_vertical_gradient_level(1) == -9999999.9 )  THEN
753          vg_vertical_gradient_level(1) = 0.0
754       ENDIF
755
756!
757!--    Compute initial temperature profile using the given temperature gradients
758       i = 1
759       gradient = 0.0
760
761       IF ( .NOT. ocean )  THEN
762
763          pt_vertical_gradient_level_ind(1) = 0
764          DO  k = 1, nzt+1
765             IF ( pt_vertical_gradient_level(i) < zu(k)  .AND. &
766                  pt_vertical_gradient_level(i) >= 0.0 )  THEN
767                gradient = pt_vertical_gradient(i) / 100.0
768                pt_vertical_gradient_level_ind(i) = k - 1
769                i = i + 1
770                IF ( i > 10 )  THEN
771                   IF ( myid == 0 )  THEN
772                      PRINT*, '+++ check_parameters: upper bound 10 of array', &
773                              ' "pt_vertical_gradient_level_ind" exceeded'
774                   ENDIF
775                   CALL local_stop
776                ENDIF
777             ENDIF
778             IF ( gradient /= 0.0 )  THEN
779                IF ( k /= 1 )  THEN
780                   pt_init(k) = pt_init(k-1) + dzu(k) * gradient
781                ELSE
782                   pt_init(k) = pt_surface   + 0.5 * dzu(k) * gradient
783                ENDIF
784             ELSE
785                pt_init(k) = pt_init(k-1)
786             ENDIF
787          ENDDO
788
789       ELSE
790
791          pt_vertical_gradient_level_ind(1) = nzt+1
792          DO  k = nzt, 0, -1
793             IF ( pt_vertical_gradient_level(i) > zu(k)  .AND. &
794                  pt_vertical_gradient_level(i) <= 0.0 )  THEN
795                gradient = pt_vertical_gradient(i) / 100.0
796                pt_vertical_gradient_level_ind(i) = k + 1
797                i = i + 1
798                IF ( i > 10 )  THEN
799                   IF ( myid == 0 )  THEN
800                      PRINT*, '+++ check_parameters: upper bound 10 of array', &
801                              ' "pt_vertical_gradient_level_ind" exceeded'
802                   ENDIF
803                   CALL local_stop
804                ENDIF
805             ENDIF
806             IF ( gradient /= 0.0 )  THEN
807                IF ( k /= nzt )  THEN
808                   pt_init(k) = pt_init(k+1) - dzu(k+1) * gradient
809                ELSE
810                   pt_init(k)   = pt_surface - 0.5 * dzu(k+1) * gradient
811                   pt_init(k+1) = pt_surface + 0.5 * dzu(k+1) * gradient
812                ENDIF
813             ELSE
814                pt_init(k) = pt_init(k+1)
815             ENDIF
816          ENDDO
817
818       ENDIF
819
820!
821!--    In case of no given temperature gradients, choose gradient of neutral
822!--    stratification
823       IF ( pt_vertical_gradient_level(1) == -9999999.9 )  THEN
824          pt_vertical_gradient_level(1) = 0.0
825       ENDIF
826
827!
828!--    Store temperature gradient at the top boundary for possible Neumann
829!--    boundary condition
830       bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1)
831
832!
833!--    If required, compute initial humidity or scalar profile using the given
834!--    humidity/scalar gradient. In case of scalar transport, initially store
835!--    values of the scalar parameters on humidity parameters
836       IF ( passive_scalar )  THEN
837          bc_q_b                    = bc_s_b
838          bc_q_t                    = bc_s_t
839          q_surface                 = s_surface
840          q_surface_initial_change  = s_surface_initial_change
841          q_vertical_gradient       = s_vertical_gradient
842          q_vertical_gradient_level = s_vertical_gradient_level
843          surface_waterflux         = surface_scalarflux
844       ENDIF
845
846       IF ( humidity  .OR.  passive_scalar )  THEN
847
848          i = 1
849          gradient = 0.0
850          q_vertical_gradient_level_ind(1) = 0
851          DO  k = 1, nzt+1
852             IF ( q_vertical_gradient_level(i) < zu(k)  .AND. &
853                  q_vertical_gradient_level(i) >= 0.0 )  THEN
854                gradient = q_vertical_gradient(i) / 100.0
855                q_vertical_gradient_level_ind(i) = k - 1
856                i = i + 1
857                IF ( i > 10 )  THEN
858                   IF ( myid == 0 )  THEN
859                      PRINT*, '+++ check_parameters: upper bound 10 of arr', &
860                              'ay "q_vertical_gradient_level_ind" exceeded'
861                   ENDIF
862                   CALL local_stop
863                ENDIF
864             ENDIF
865             IF ( gradient /= 0.0 )  THEN
866                IF ( k /= 1 )  THEN
867                   q_init(k) = q_init(k-1) + dzu(k) * gradient
868                ELSE
869                   q_init(k) = q_init(k-1) + 0.5 * dzu(k) * gradient
870                ENDIF
871             ELSE
872                q_init(k) = q_init(k-1)
873             ENDIF
874!
875!--          Avoid negative humidities
876             IF ( q_init(k) < 0.0 )  THEN
877                q_init(k) = 0.0
878             ENDIF
879          ENDDO
880
881!
882!--       In case of no given humidity gradients, choose zero gradient
883!--       conditions
884          IF ( q_vertical_gradient_level(1) == -1.0 )  THEN
885             q_vertical_gradient_level(1) = 0.0
886          ENDIF
887
888!
889!--       Store humidity gradient at the top boundary for possile Neumann
890!--       boundary condition
891          bc_q_t_val = ( q_init(nzt+1) - q_init(nzt) ) / dzu(nzt+1)
892
893       ENDIF
894
895!
896!--    If required, compute initial salinity profile using the given salinity
897!--    gradients
898       IF ( ocean )  THEN
899
900          i = 1
901          gradient = 0.0
902
903          sa_vertical_gradient_level_ind(1) = nzt+1
904          DO  k = nzt, 0, -1
905             IF ( sa_vertical_gradient_level(i) > zu(k)  .AND. &
906                  sa_vertical_gradient_level(i) <= 0.0 )  THEN
907                gradient = sa_vertical_gradient(i) / 100.0
908                sa_vertical_gradient_level_ind(i) = k + 1
909                i = i + 1
910                IF ( i > 10 )  THEN
911                   IF ( myid == 0 )  THEN
912                      PRINT*, '+++ check_parameters: upper bound 10 of array', &
913                              ' "sa_vertical_gradient_level_ind" exceeded'
914                   ENDIF
915                   CALL local_stop
916                ENDIF
917             ENDIF
918             IF ( gradient /= 0.0 )  THEN
919                IF ( k /= nzt )  THEN
920                   sa_init(k) = sa_init(k+1) - dzu(k+1) * gradient
921                ELSE
922                   sa_init(k)   = sa_surface - 0.5 * dzu(k+1) * gradient
923                   sa_init(k+1) = sa_surface + 0.5 * dzu(k+1) * gradient
924                ENDIF
925             ELSE
926                sa_init(k) = sa_init(k+1)
927             ENDIF
928          ENDDO
929
930       ENDIF
931
932    ENDIF
933
934!
935!-- Compute Coriolis parameter
936    f  = 2.0 * omega * SIN( phi / 180.0 * pi )
937    fs = 2.0 * omega * COS( phi / 180.0 * pi )
938
939!
940!-- Ocean runs always use reference values in the buoyancy term. Therefore
941!-- set the reference temperature equal to the surface temperature.
942    IF ( ocean  .AND.  pt_reference == 9999999.9 )  pt_reference = pt_surface
943
944!
945!-- Reference value has to be used in buoyancy terms
946    IF ( pt_reference /= 9999999.9 )  use_reference = .TRUE.
947
948!
949!-- Sign of buoyancy/stability terms
950    IF ( ocean )  atmos_ocean_sign = -1.0
951
952!
953!-- Ocean version must use flux boundary conditions at the top
954    IF ( ocean .AND. .NOT. use_top_fluxes )  THEN
955       IF ( myid == 0 )  PRINT*, '+++ check_parameters: use_top_fluxes ',&
956                                    'must be .TRUE. in ocean version'
957       CALL local_stop
958    ENDIF
959
960!
961!-- In case of a given slope, compute the relevant quantities
962    IF ( alpha_surface /= 0.0 )  THEN
963       IF ( ABS( alpha_surface ) > 90.0 )  THEN
964          IF ( myid == 0 )  PRINT*, '+++ check_parameters: ABS( alpha_surface',&
965                                    '=', alpha_surface, ' ) must be < 90.0'
966          CALL local_stop
967       ENDIF
968       sloping_surface = .TRUE.
969       cos_alpha_surface = COS( alpha_surface / 180.0 * pi )
970       sin_alpha_surface = SIN( alpha_surface / 180.0 * pi )
971    ENDIF
972
973!
974!-- Check time step and cfl_factor
975    IF ( dt /= -1.0 )  THEN
976       IF ( dt <= 0.0  .AND.  dt /= -1.0 )  THEN
977          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  dt=', dt, ' <= 0.0'
978          CALL local_stop
979       ENDIF
980       dt_3d = dt
981       dt_fixed = .TRUE.
982    ENDIF
983
984    IF ( cfl_factor <= 0.0  .OR.  cfl_factor > 1.0 )  THEN
985       IF ( cfl_factor == -1.0 )  THEN
986          IF ( momentum_advec == 'ups-scheme'  .OR.  &
987               scalar_advec == 'ups-scheme' )  THEN
988             cfl_factor = 0.8
989          ELSE
990             IF ( timestep_scheme == 'runge-kutta-2' )  THEN
991                cfl_factor = 0.8
992             ELSEIF ( timestep_scheme == 'runge-kutta-3' )  THEN
993                cfl_factor = 0.9
994             ELSE
995                cfl_factor = 0.1
996             ENDIF
997          ENDIF
998       ELSE
999          IF ( myid == 0 )  THEN
1000             PRINT*, '+++ check_parameters: cfl_factor=', cfl_factor, &
1001                         ' out of range'
1002             PRINT*, '+++                   0.0 < cfl_factor <= 1.0 is required'
1003          ENDIF
1004          CALL local_stop
1005       ENDIF
1006    ENDIF
1007
1008!
1009!-- Store simulated time at begin
1010    simulated_time_at_begin = simulated_time
1011
1012!
1013!-- Set wind speed in the Galilei-transformed system
1014    IF ( galilei_transformation )  THEN
1015       IF ( use_ug_for_galilei_tr .AND.                &
1016            ug_vertical_gradient_level(1) == 0.0 .AND. & 
1017            vg_vertical_gradient_level(1) == 0.0 )  THEN
1018          u_gtrans = ug_surface
1019          v_gtrans = vg_surface
1020       ELSEIF ( use_ug_for_galilei_tr .AND.                &
1021                ug_vertical_gradient_level(1) /= 0.0 )  THEN
1022          IF ( myid == 0 )  THEN
1023             PRINT*, '+++ check_parameters:'
1024             PRINT*, '    baroclinicity (ug) not allowed'
1025             PRINT*, '    simultaneously with galilei transformation'
1026          ENDIF
1027          CALL local_stop
1028       ELSEIF ( use_ug_for_galilei_tr .AND.                &
1029                vg_vertical_gradient_level(1) /= 0.0 )  THEN
1030          IF ( myid == 0 )  THEN
1031             PRINT*, '+++ check_parameters:'
1032             PRINT*, '    baroclinicity (vg) not allowed'
1033             PRINT*, '    simultaneously with galilei transformation'
1034          ENDIF
1035          CALL local_stop
1036       ELSE
1037          IF ( myid == 0 )  THEN
1038             PRINT*, '+++ WARNING: check_parameters:'
1039             PRINT*, '    variable translation speed used for galilei-tran' // &
1040                          'sformation, which'
1041             PRINT*, '    may cause instabilities in stably stratified regions'
1042          ENDIF
1043       ENDIF
1044    ENDIF
1045
1046!
1047!-- In case of using a prandtl-layer, calculated (or prescribed) surface
1048!-- fluxes have to be used in the diffusion-terms
1049    IF ( prandtl_layer )  use_surface_fluxes = .TRUE.
1050
1051!
1052!-- Check boundary conditions and set internal variables:
1053!-- Lateral boundary conditions
1054    IF ( bc_lr /= 'cyclic'  .AND.  bc_lr /= 'dirichlet/radiation'  .AND. &
1055         bc_lr /= 'radiation/dirichlet' )  THEN
1056       IF ( myid == 0 )  THEN
1057          PRINT*, '+++ check_parameters:'
1058          PRINT*, '    unknown boundary condition: bc_lr = ', bc_lr
1059       ENDIF
1060       CALL local_stop
1061    ENDIF
1062    IF ( bc_ns /= 'cyclic'  .AND.  bc_ns /= 'dirichlet/radiation'  .AND. &
1063         bc_ns /= 'radiation/dirichlet' )  THEN
1064       IF ( myid == 0 )  THEN
1065          PRINT*, '+++ check_parameters:'
1066          PRINT*, '    unknown boundary condition: bc_ns = ', bc_ns
1067       ENDIF
1068       CALL local_stop
1069    ENDIF
1070
1071!
1072!-- Non-cyclic lateral boundaries require the multigrid method and Piascek-
1073!-- Willimas advection scheme. Several schemes and tools do not work with
1074!-- non-cyclic boundary conditions.
1075    IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
1076       IF ( psolver /= 'multigrid' )  THEN
1077          IF ( myid == 0 )  THEN
1078             PRINT*, '+++ check_parameters:'
1079             PRINT*, '    non-cyclic lateral boundaries do not allow', &
1080                          ' psolver = ', psolver 
1081          ENDIF
1082          CALL local_stop
1083       ENDIF
1084       IF ( momentum_advec /= 'pw-scheme' )  THEN
1085          IF ( myid == 0 )  THEN
1086             PRINT*, '+++ check_parameters:'
1087             PRINT*, '    non-cyclic lateral boundaries do not allow', &
1088                          ' momentum_advec = ', momentum_advec 
1089          ENDIF
1090          CALL local_stop
1091       ENDIF
1092       IF ( scalar_advec /= 'pw-scheme' )  THEN
1093          IF ( myid == 0 )  THEN
1094             PRINT*, '+++ check_parameters:'
1095             PRINT*, '    non-cyclic lateral boundaries do not allow', &
1096                          ' scalar_advec = ', scalar_advec 
1097          ENDIF
1098          CALL local_stop
1099       ENDIF
1100       IF ( galilei_transformation )  THEN
1101          IF ( myid == 0 )  THEN
1102             PRINT*, '+++ check_parameters:'
1103             PRINT*, '    non-cyclic lateral boundaries do not allow', &
1104                          ' galilei_transformation = .T.' 
1105          ENDIF
1106          CALL local_stop
1107       ENDIF
1108!       IF ( conserve_volume_flow )  THEN
1109!          IF ( myid == 0 )  THEN
1110!             PRINT*, '+++ check_parameters:'
1111!             PRINT*, '    non-cyclic lateral boundaries do not allow', &
1112!                          ' conserve_volume_flow = .T.'
1113!          ENDIF
1114!          CALL local_stop
1115!       ENDIF
1116    ENDIF
1117
1118!
1119!-- Bottom boundary condition for the turbulent Kinetic energy
1120    IF ( bc_e_b == 'neumann' )  THEN
1121       ibc_e_b = 1
1122       IF ( adjust_mixing_length  .AND.  prandtl_layer )  THEN
1123          IF ( myid == 0 )  THEN
1124             PRINT*, '+++ WARNING: check_parameters:'
1125             PRINT*, '    adjust_mixing_length = TRUE and bc_e_b = ', bc_e_b
1126          ENDIF
1127       ENDIF
1128    ELSEIF ( bc_e_b == '(u*)**2+neumann' )  THEN
1129       ibc_e_b = 2
1130       IF ( .NOT. adjust_mixing_length  .AND.  prandtl_layer )  THEN
1131          IF ( myid == 0 )  THEN
1132             PRINT*, '+++ WARNING: check_parameters:'
1133             PRINT*, '    adjust_mixing_length = FALSE and bc_e_b = ', bc_e_b
1134          ENDIF
1135       ENDIF
1136       IF ( .NOT. prandtl_layer )  THEN
1137          bc_e_b = 'neumann'
1138          ibc_e_b = 1
1139          IF ( myid == 0 )  THEN
1140             PRINT*, '+++ WARNING: check_parameters:'
1141             PRINT*, '    boundary condition bc_e_b changed to "', bc_e_b, '"'
1142          ENDIF
1143       ENDIF
1144    ELSE
1145       IF ( myid == 0 )  THEN
1146          PRINT*, '+++ check_parameters:'
1147          PRINT*, '    unknown boundary condition: bc_e_b = ', bc_e_b
1148       ENDIF
1149       CALL local_stop
1150    ENDIF
1151
1152!
1153!-- Boundary conditions for perturbation pressure
1154    IF ( bc_p_b == 'dirichlet' )  THEN
1155       ibc_p_b = 0
1156    ELSEIF ( bc_p_b == 'neumann' )  THEN
1157       ibc_p_b = 1
1158    ELSEIF ( bc_p_b == 'neumann+inhomo' )  THEN
1159       ibc_p_b = 2
1160    ELSE
1161       IF ( myid == 0 )  THEN
1162          PRINT*, '+++ check_parameters:'
1163          PRINT*, '    unknown boundary condition: bc_p_b = ', bc_p_b
1164       ENDIF
1165       CALL local_stop
1166    ENDIF
1167    IF ( ibc_p_b == 2  .AND.  .NOT. prandtl_layer )  THEN
1168       IF ( myid == 0 )  THEN
1169          PRINT*, '+++ check_parameters:'
1170          PRINT*, '    boundary condition: bc_p_b = ', TRIM( bc_p_b ), &
1171                       ' not allowed with'
1172          PRINT*, '    prandtl_layer = .FALSE.' 
1173       ENDIF
1174       CALL local_stop
1175    ENDIF
1176    IF ( bc_p_t == 'dirichlet' )  THEN
1177       ibc_p_t = 0
1178    ELSEIF ( bc_p_t == 'neumann' )  THEN
1179       ibc_p_t = 1
1180    ELSE
1181       IF ( myid == 0 )  THEN
1182          PRINT*, '+++ check_parameters:'
1183          PRINT*, '    unknown boundary condition: bc_p_t = ', bc_p_t
1184       ENDIF
1185       CALL local_stop
1186    ENDIF
1187
1188!
1189!-- Boundary conditions for potential temperature
1190    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
1191       ibc_pt_b = 2
1192    ELSE
1193       IF ( bc_pt_b == 'dirichlet' )  THEN
1194          ibc_pt_b = 0
1195       ELSEIF ( bc_pt_b == 'neumann' )  THEN
1196          ibc_pt_b = 1
1197       ELSE
1198          IF ( myid == 0 )  THEN
1199             PRINT*, '+++ check_parameters:'
1200             PRINT*, '    unknown boundary condition: bc_pt_b = ', bc_pt_b
1201          ENDIF
1202          CALL local_stop
1203       ENDIF
1204    ENDIF
1205
1206    IF ( bc_pt_t == 'dirichlet' )  THEN
1207       ibc_pt_t = 0
1208    ELSEIF ( bc_pt_t == 'neumann' )  THEN
1209       ibc_pt_t = 1
1210    ELSEIF ( bc_pt_t == 'initial_gradient' )  THEN
1211       ibc_pt_t = 2
1212    ELSE
1213       IF ( myid == 0 )  THEN
1214          PRINT*, '+++ check_parameters:'
1215          PRINT*, '    unknown boundary condition: bc_pt_t = ', bc_pt_t
1216       ENDIF
1217       CALL local_stop
1218    ENDIF
1219
1220    IF ( surface_heatflux == 9999999.9 )  constant_heatflux     = .FALSE.
1221    IF ( top_heatflux     == 9999999.9 )  constant_top_heatflux = .FALSE.
1222    IF ( top_momentumflux_u /= 9999999.9  .AND.  &
1223         top_momentumflux_v /= 9999999.9 )  THEN
1224       constant_top_momentumflux = .TRUE.
1225    ELSEIF (  .NOT. ( top_momentumflux_u == 9999999.9  .AND.  &
1226           top_momentumflux_v == 9999999.9 ) )  THEN   
1227       IF ( myid == 0 )  THEN
1228          PRINT*, '+++ check_parameters:'
1229          PRINT*, '    both, top_momentumflux_u AND top_momentumflux_v'
1230          PRINT*, '    must be set'
1231       ENDIF
1232       CALL local_stop
1233    ENDIF
1234
1235!
1236!-- A given surface temperature implies Dirichlet boundary condition for
1237!-- temperature. In this case specification of a constant heat flux is
1238!-- forbidden.
1239    IF ( ibc_pt_b == 0  .AND.   constant_heatflux  .AND. &
1240         surface_heatflux /= 0.0 )  THEN
1241       IF ( myid == 0 )  THEN
1242          PRINT*, '+++ check_parameters:'
1243          PRINT*, '    boundary_condition: bc_pt_b = ', bc_pt_b
1244          PRINT*, '    is not allowed with constant_heatflux = .TRUE.'
1245       ENDIF
1246       CALL local_stop
1247    ENDIF
1248    IF ( constant_heatflux  .AND.  pt_surface_initial_change /= 0.0 )  THEN
1249       IF ( myid == 0 )  THEN
1250          PRINT*, '+++ check_parameters: constant_heatflux = .TRUE. is not'
1251          PRINT*, '    allowed with pt_surface_initial_change (/=0) = ', &
1252                  pt_surface_initial_change
1253       ENDIF
1254       CALL local_stop
1255    ENDIF
1256
1257!
1258!-- A given temperature at the top implies Dirichlet boundary condition for
1259!-- temperature. In this case specification of a constant heat flux is
1260!-- forbidden.
1261    IF ( ibc_pt_t == 0  .AND.   constant_top_heatflux  .AND. &
1262         top_heatflux /= 0.0 )  THEN
1263       IF ( myid == 0 )  THEN
1264          PRINT*, '+++ check_parameters:'
1265          PRINT*, '    boundary_condition: bc_pt_t = ', bc_pt_t
1266          PRINT*, '    is not allowed with constant_top_heatflux = .TRUE.'
1267       ENDIF
1268       CALL local_stop
1269    ENDIF
1270
1271!
1272!-- Boundary conditions for salinity
1273    IF ( ocean )  THEN
1274       IF ( bc_sa_t == 'dirichlet' )  THEN
1275          ibc_sa_t = 0
1276       ELSEIF ( bc_sa_t == 'neumann' )  THEN
1277          ibc_sa_t = 1
1278       ELSE
1279          IF ( myid == 0 )  THEN
1280             PRINT*, '+++ check_parameters:'
1281             PRINT*, '    unknown boundary condition: bc_sa_t = ', bc_sa_t
1282          ENDIF
1283          CALL local_stop
1284       ENDIF
1285
1286       IF ( top_salinityflux == 9999999.9 )  constant_top_salinityflux = .FALSE.
1287       IF ( ibc_sa_t == 1  .AND.   top_salinityflux == 9999999.9 )  THEN
1288          IF ( myid == 0 )  THEN
1289             PRINT*, '+++ check_parameters:'
1290             PRINT*, '    boundary_condition: bc_sa_t = ', bc_sa_t
1291             PRINT*, '    requires to set top_salinityflux '
1292          ENDIF
1293          CALL local_stop
1294       ENDIF
1295
1296!
1297!--    A fixed salinity at the top implies Dirichlet boundary condition for
1298!--    salinity. In this case specification of a constant salinity flux is
1299!--    forbidden.
1300       IF ( ibc_sa_t == 0  .AND.   constant_top_salinityflux  .AND. &
1301            top_salinityflux /= 0.0 )  THEN
1302          IF ( myid == 0 )  THEN
1303             PRINT*, '+++ check_parameters:'
1304             PRINT*, '    boundary_condition: bc_sa_t = ', bc_sa_t
1305             PRINT*, '    is not allowed with constant_top_salinityflux = ', &
1306                          '.TRUE.'
1307          ENDIF
1308          CALL local_stop
1309       ENDIF
1310
1311    ENDIF
1312
1313!
1314!-- In case of humidity or passive scalar, set boundary conditions for total
1315!-- water content / scalar
1316    IF ( humidity  .OR.  passive_scalar ) THEN
1317       IF ( humidity )  THEN
1318          sq = 'q'
1319       ELSE
1320          sq = 's'
1321       ENDIF
1322       IF ( bc_q_b == 'dirichlet' )  THEN
1323          ibc_q_b = 0
1324       ELSEIF ( bc_q_b == 'neumann' )  THEN
1325          ibc_q_b = 1
1326       ELSE
1327          IF ( myid == 0 )  THEN
1328             PRINT*, '+++ check_parameters:'
1329             PRINT*, '    unknown boundary condition: bc_', sq, '_b = ', bc_q_b
1330          ENDIF
1331          CALL local_stop
1332       ENDIF
1333       IF ( bc_q_t == 'dirichlet' )  THEN
1334          ibc_q_t = 0
1335       ELSEIF ( bc_q_t == 'neumann' )  THEN
1336          ibc_q_t = 1
1337       ELSE
1338          IF ( myid == 0 )  THEN
1339             PRINT*, '+++ check_parameters:'
1340             PRINT*, '    unknown boundary condition: bc_', sq, '_t = ', bc_q_t
1341          ENDIF
1342          CALL local_stop
1343       ENDIF
1344
1345       IF ( surface_waterflux == 0.0 )  constant_waterflux = .FALSE.
1346
1347!
1348!--    A given surface humidity implies Dirichlet boundary condition for
1349!--    humidity. In this case specification of a constant water flux is
1350!--    forbidden.
1351       IF ( ibc_q_b == 0  .AND.  constant_waterflux )  THEN
1352          IF ( myid == 0 )  THEN
1353             PRINT*, '+++ check_parameters:'
1354             PRINT*, '    boundary_condition: bc_', sq, '_b = ', bc_q_b
1355             PRINT*, '    is not allowed with prescribed surface flux'
1356          ENDIF
1357          CALL local_stop
1358       ENDIF
1359       IF ( constant_waterflux  .AND.  q_surface_initial_change /= 0.0 )  THEN
1360          IF ( myid == 0 )  THEN
1361             PRINT*, '+++ check_parameters: a prescribed surface flux is not'
1362             PRINT*, '    allowed with ', sq, '_surface_initial_change (/=0)', &
1363                     ' = ', q_surface_initial_change
1364          ENDIF
1365          CALL local_stop
1366       ENDIF
1367       
1368    ENDIF
1369
1370!
1371!-- Boundary conditions for horizontal components of wind speed
1372    IF ( bc_uv_b == 'dirichlet' )  THEN
1373       ibc_uv_b = 0
1374    ELSEIF ( bc_uv_b == 'neumann' )  THEN
1375       ibc_uv_b = 1
1376       IF ( prandtl_layer )  THEN
1377          IF ( myid == 0 )  THEN
1378             PRINT*, '+++ check_parameters:'
1379             PRINT*, '    boundary condition: bc_uv_b = ', TRIM( bc_uv_b ), &
1380                          ' is not allowed with'
1381             PRINT*, '    prandtl_layer = .TRUE.' 
1382          ENDIF
1383          CALL local_stop
1384       ENDIF
1385    ELSE
1386       IF ( myid == 0 )  THEN
1387          PRINT*, '+++ check_parameters:'
1388          PRINT*, '    unknown boundary condition: bc_uv_b = ', bc_uv_b
1389       ENDIF
1390       CALL local_stop
1391    ENDIF
1392    IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1393       bc_uv_t = 'neumann'
1394       ibc_uv_t = 1
1395    ELSE
1396       IF ( bc_uv_t == 'dirichlet' )  THEN
1397          ibc_uv_t = 0
1398       ELSEIF ( bc_uv_t == 'neumann' )  THEN
1399          ibc_uv_t = 1
1400       ELSE
1401          IF ( myid == 0 )  THEN
1402             PRINT*, '+++ check_parameters:'
1403             PRINT*, '    unknown boundary condition: bc_uv_t = ', bc_uv_t
1404          ENDIF
1405          CALL local_stop
1406       ENDIF
1407    ENDIF
1408
1409!
1410!-- Compute and check, respectively, the Rayleigh Damping parameter
1411    IF ( rayleigh_damping_factor == -1.0 )  THEN
1412       IF ( momentum_advec == 'ups-scheme' )  THEN
1413          rayleigh_damping_factor = 0.01
1414       ELSE
1415          rayleigh_damping_factor = 0.0
1416       ENDIF
1417    ELSE
1418       IF ( rayleigh_damping_factor < 0.0 .OR. rayleigh_damping_factor > 1.0 ) &
1419       THEN
1420          IF ( myid == 0 )  THEN
1421             PRINT*, '+++ check_parameters:'
1422             PRINT*, '    rayleigh_damping_factor = ', rayleigh_damping_factor,&
1423                          ' out of range [0.0,1.0]'
1424          ENDIF
1425          CALL local_stop
1426       ENDIF
1427    ENDIF
1428
1429    IF ( rayleigh_damping_height == -1.0 )  THEN
1430       IF ( .NOT. ocean )  THEN
1431          rayleigh_damping_height = 0.66666666666 * zu(nzt)
1432       ELSE
1433          rayleigh_damping_height = 0.66666666666 * zu(nzb)
1434       ENDIF
1435    ELSE
1436       IF ( .NOT. ocean )  THEN
1437          IF ( rayleigh_damping_height < 0.0  .OR. &
1438               rayleigh_damping_height > zu(nzt) )  THEN
1439             IF ( myid == 0 )  THEN
1440                PRINT*, '+++ check_parameters:'
1441                PRINT*, '    rayleigh_damping_height = ', rayleigh_damping_height,&
1442                     ' out of range [0.0,', zu(nzt), ']'
1443             ENDIF
1444             CALL local_stop
1445          ENDIF
1446       ELSE
1447          IF ( rayleigh_damping_height > 0.0  .OR. &
1448               rayleigh_damping_height < zu(nzb) )  THEN
1449             IF ( myid == 0 )  THEN
1450                PRINT*, '+++ check_parameters:'
1451                PRINT*, '    rayleigh_damping_height = ', rayleigh_damping_height,&
1452                     ' out of range [0.0,', zu(nzb), ']'
1453             ENDIF
1454             CALL local_stop
1455          ENDIF
1456       ENDIF
1457    ENDIF
1458
1459!
1460!-- Check limiters for Upstream-Spline scheme
1461    IF ( overshoot_limit_u < 0.0  .OR.  overshoot_limit_v < 0.0  .OR.  &
1462         overshoot_limit_w < 0.0  .OR.  overshoot_limit_pt < 0.0  .OR. &
1463         overshoot_limit_e < 0.0 )  THEN
1464       IF ( myid == 0 )  THEN
1465          PRINT*, '+++ check_parameters:'
1466          PRINT*, '    overshoot_limit_... < 0.0 is not allowed'
1467       ENDIF
1468       CALL local_stop
1469    ENDIF
1470    IF ( ups_limit_u < 0.0 .OR. ups_limit_v < 0.0 .OR. ups_limit_w < 0.0 .OR. &
1471         ups_limit_pt < 0.0 .OR. ups_limit_e < 0.0 )  THEN
1472       IF ( myid == 0 )  THEN
1473          PRINT*, '+++ check_parameters:'
1474          PRINT*, '    ups_limit_... < 0.0 is not allowed'
1475       ENDIF
1476       CALL local_stop
1477    ENDIF
1478
1479!
1480!-- Check number of chosen statistic regions. More than 10 regions are not
1481!-- allowed, because so far no more than 10 corresponding output files can
1482!-- be opened (cf. check_open)
1483    IF ( statistic_regions > 9  .OR.  statistic_regions < 0 )  THEN
1484       IF ( myid == 0 )  THEN
1485          PRINT*, '+++ check_parameters: Number of statistic_regions = ', &
1486                       statistic_regions+1
1487          PRINT*, '    Only 10 regions are allowed'
1488       ENDIF
1489       CALL local_stop
1490    ENDIF
1491    IF ( normalizing_region > statistic_regions  .OR. &
1492         normalizing_region < 0)  THEN
1493       IF ( myid == 0 )  THEN
1494          PRINT*, '+++ check_parameters: normalizing_region = ', &
1495                       normalizing_region, ' is unknown'
1496          PRINT*, '    Must be <= ', statistic_regions
1497       ENDIF
1498       CALL local_stop
1499    ENDIF
1500
1501!
1502!-- Set the default intervals for data output, if necessary
1503!-- NOTE: dt_dosp has already been set in package_parin
1504    IF ( dt_data_output /= 9999999.9 )  THEN
1505       IF ( dt_dopr           == 9999999.9 )  dt_dopr           = dt_data_output
1506       IF ( dt_dopts          == 9999999.9 )  dt_dopts          = dt_data_output
1507       IF ( dt_do2d_xy        == 9999999.9 )  dt_do2d_xy        = dt_data_output
1508       IF ( dt_do2d_xz        == 9999999.9 )  dt_do2d_xz        = dt_data_output
1509       IF ( dt_do2d_yz        == 9999999.9 )  dt_do2d_yz        = dt_data_output
1510       IF ( dt_do3d           == 9999999.9 )  dt_do3d           = dt_data_output
1511       IF ( dt_data_output_av == 9999999.9 )  dt_data_output_av = dt_data_output
1512    ENDIF
1513
1514!
1515!-- Set the default skip time intervals for data output, if necessary
1516    IF ( skip_time_dopr    == 9999999.9 ) &
1517                                       skip_time_dopr    = skip_time_data_output
1518    IF ( skip_time_dosp    == 9999999.9 ) &
1519                                       skip_time_dosp    = skip_time_data_output
1520    IF ( skip_time_do2d_xy == 9999999.9 ) &
1521                                       skip_time_do2d_xy = skip_time_data_output
1522    IF ( skip_time_do2d_xz == 9999999.9 ) &
1523                                       skip_time_do2d_xz = skip_time_data_output
1524    IF ( skip_time_do2d_yz == 9999999.9 ) &
1525                                       skip_time_do2d_yz = skip_time_data_output
1526    IF ( skip_time_do3d    == 9999999.9 ) &
1527                                       skip_time_do3d    = skip_time_data_output
1528    IF ( skip_time_data_output_av == 9999999.9 ) &
1529                                skip_time_data_output_av = skip_time_data_output
1530
1531!
1532!-- Check the average intervals (first for 3d-data, then for profiles and
1533!-- spectra)
1534    IF ( averaging_interval > dt_data_output_av )  THEN
1535       IF ( myid == 0 )  THEN
1536          PRINT*, '+++ check_parameters: average_interval=',              &
1537                       averaging_interval, ' must be <= dt_data_output=', &
1538                       dt_data_output
1539       ENDIF
1540       CALL local_stop
1541    ENDIF
1542
1543    IF ( averaging_interval_pr == 9999999.9 )  THEN
1544       averaging_interval_pr = averaging_interval
1545    ENDIF
1546
1547    IF ( averaging_interval_pr > dt_dopr )  THEN
1548       IF ( myid == 0 )  THEN
1549          PRINT*, '+++ check_parameters: averaging_interval_pr=', &
1550                       averaging_interval_pr, ' must be <= dt_dopr=', dt_dopr
1551       ENDIF
1552       CALL local_stop
1553    ENDIF
1554
1555    IF ( averaging_interval_sp == 9999999.9 )  THEN
1556       averaging_interval_sp = averaging_interval
1557    ENDIF
1558
1559    IF ( averaging_interval_sp > dt_dosp )  THEN
1560       IF ( myid == 0 )  THEN
1561          PRINT*, '+++ check_parameters: averaging_interval_sp=', &
1562                       averaging_interval_sp, ' must be <= dt_dosp=', dt_dosp
1563       ENDIF
1564       CALL local_stop
1565    ENDIF
1566
1567!
1568!-- Set the default interval for profiles entering the temporal average
1569    IF ( dt_averaging_input_pr == 9999999.9 )  THEN
1570       dt_averaging_input_pr = dt_averaging_input
1571    ENDIF
1572
1573!
1574!-- Set the default interval for the output of timeseries to a reasonable
1575!-- value (tries to minimize the number of calls of flow_statistics)
1576    IF ( dt_dots == 9999999.9 )  THEN
1577       IF ( averaging_interval_pr == 0.0 )  THEN
1578          dt_dots = MIN( dt_run_control, dt_dopr )
1579       ELSE
1580          dt_dots = MIN( dt_run_control, dt_averaging_input_pr )
1581       ENDIF
1582    ENDIF
1583
1584!
1585!-- Check the sample rate for averaging (first for 3d-data, then for profiles)
1586    IF ( dt_averaging_input > averaging_interval )  THEN
1587       IF ( myid == 0 )  THEN
1588          PRINT*, '+++ check_parameters: dt_averaging_input=',                &
1589                       dt_averaging_input, ' must be <= averaging_interval=', &
1590                       averaging_interval
1591       ENDIF
1592       CALL local_stop
1593    ENDIF
1594
1595    IF ( dt_averaging_input_pr > averaging_interval_pr )  THEN
1596       IF ( myid == 0 )  THEN
1597          PRINT*, '+++ check_parameters: dt_averaging_input_pr=', &
1598                       dt_averaging_input_pr,                     &
1599                       ' must be <= averaging_interval_pr=',      &
1600                       averaging_interval_pr
1601       ENDIF
1602       CALL local_stop
1603    ENDIF
1604
1605!
1606!-- Set the default value for the integration interval of precipitation amount
1607    IF ( precipitation )  THEN
1608       IF ( precipitation_amount_interval == 9999999.9 )  THEN
1609          precipitation_amount_interval = dt_do2d_xy
1610       ELSE
1611          IF ( precipitation_amount_interval > dt_do2d_xy )  THEN
1612             IF ( myid == 0 )  PRINT*, '+++ check_parameters: ',              &
1613                                       'precipitation_amount_interval =',     &
1614                                        precipitation_amount_interval,        &
1615                                       ' must not be larger than dt_do2d_xy', &
1616                                       ' = ', dt_do2d_xy   
1617       CALL local_stop
1618          ENDIF
1619       ENDIF
1620    ENDIF
1621
1622!
1623!-- Determine the number of output profiles and check whether they are
1624!-- permissible
1625    DO  WHILE ( data_output_pr(dopr_n+1) /= '          ' )
1626
1627       dopr_n = dopr_n + 1
1628       i = dopr_n
1629
1630!
1631!--    Determine internal profile number (for hom, homs)
1632!--    and store height levels
1633       SELECT CASE ( TRIM( data_output_pr(i) ) )
1634
1635          CASE ( 'u', '#u' )
1636             dopr_index(i) = 1
1637             dopr_unit(i)  = 'm/s'
1638             hom(:,2,1,:)  = SPREAD( zu, 2, statistic_regions+1 )
1639             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1640                dopr_initial_index(i) = 5
1641                hom(:,2,5,:)          = SPREAD( zu, 2, statistic_regions+1 )
1642                data_output_pr(i)     = data_output_pr(i)(2:)
1643             ENDIF
1644
1645          CASE ( 'v', '#v' )
1646             dopr_index(i) = 2
1647             dopr_unit(i)  = 'm/s'
1648             hom(:,2,2,:)  = SPREAD( zu, 2, statistic_regions+1 )
1649             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1650                dopr_initial_index(i) = 6
1651                hom(:,2,6,:)          = SPREAD( zu, 2, statistic_regions+1 )
1652                data_output_pr(i)     = data_output_pr(i)(2:)
1653             ENDIF
1654
1655          CASE ( 'w' )
1656             dopr_index(i) = 3
1657             dopr_unit(i)  = 'm/s'
1658             hom(:,2,3,:)  = SPREAD( zw, 2, statistic_regions+1 )
1659
1660          CASE ( 'pt', '#pt' )
1661             IF ( .NOT. cloud_physics ) THEN
1662                dopr_index(i) = 4
1663                dopr_unit(i)  = 'K'
1664                hom(:,2,4,:)  = SPREAD( zu, 2, statistic_regions+1 )
1665                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1666                   dopr_initial_index(i) = 7
1667                   hom(:,2,7,:)          = SPREAD( zu, 2, statistic_regions+1 )
1668                   hom(nzb,2,7,:)        = 0.0    ! because zu(nzb) is negative
1669                   data_output_pr(i)     = data_output_pr(i)(2:)
1670                ENDIF
1671             ELSE
1672                dopr_index(i) = 43
1673                dopr_unit(i)  = 'K'
1674                hom(:,2,43,:)  = SPREAD( zu, 2, statistic_regions+1 )
1675                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1676                   dopr_initial_index(i) = 28
1677                   hom(:,2,28,:)         = SPREAD( zu, 2, statistic_regions+1 )
1678                   hom(nzb,2,28,:)       = 0.0    ! because zu(nzb) is negative
1679                   data_output_pr(i)     = data_output_pr(i)(2:)
1680                ENDIF
1681             ENDIF
1682
1683          CASE ( 'e' )
1684             dopr_index(i)  = 8
1685             dopr_unit(i)   = 'm2/s2'
1686             hom(:,2,8,:)   = SPREAD( zu, 2, statistic_regions+1 )
1687             hom(nzb,2,8,:) = 0.0
1688
1689          CASE ( 'km', '#km' )
1690             dopr_index(i)  = 9
1691             dopr_unit(i)   = 'm2/s'
1692             hom(:,2,9,:)   = SPREAD( zu, 2, statistic_regions+1 )
1693             hom(nzb,2,9,:) = 0.0
1694             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1695                dopr_initial_index(i) = 23
1696                hom(:,2,23,:)         = hom(:,2,9,:)
1697                data_output_pr(i)     = data_output_pr(i)(2:)
1698             ENDIF
1699
1700          CASE ( 'kh', '#kh' )
1701             dopr_index(i)   = 10
1702             dopr_unit(i)    = 'm2/s'
1703             hom(:,2,10,:)   = SPREAD( zu, 2, statistic_regions+1 )
1704             hom(nzb,2,10,:) = 0.0
1705             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1706                dopr_initial_index(i) = 24
1707                hom(:,2,24,:)         = hom(:,2,10,:)
1708                data_output_pr(i)     = data_output_pr(i)(2:)
1709             ENDIF
1710
1711          CASE ( 'l', '#l' )
1712             dopr_index(i)   = 11
1713             dopr_unit(i)    = 'm'
1714             hom(:,2,11,:)   = SPREAD( zu, 2, statistic_regions+1 )
1715             hom(nzb,2,11,:) = 0.0
1716             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1717                dopr_initial_index(i) = 25
1718                hom(:,2,25,:)         = hom(:,2,11,:)
1719                data_output_pr(i)     = data_output_pr(i)(2:)
1720             ENDIF
1721
1722          CASE ( 'w"u"' )
1723             dopr_index(i) = 12
1724             dopr_unit(i)  = 'm2/s2'
1725             hom(:,2,12,:) = SPREAD( zw, 2, statistic_regions+1 )
1726             IF ( prandtl_layer )  hom(nzb,2,12,:) = zu(1)
1727
1728          CASE ( 'w*u*' )
1729             dopr_index(i) = 13
1730             dopr_unit(i)  = 'm2/s2'
1731             hom(:,2,13,:) = SPREAD( zw, 2, statistic_regions+1 )
1732
1733          CASE ( 'w"v"' )
1734             dopr_index(i) = 14
1735             dopr_unit(i)  = 'm2/s2'
1736             hom(:,2,14,:) = SPREAD( zw, 2, statistic_regions+1 )
1737             IF ( prandtl_layer )  hom(nzb,2,14,:) = zu(1)
1738
1739          CASE ( 'w*v*' )
1740             dopr_index(i) = 15
1741             dopr_unit(i)  = 'm2/s2'
1742             hom(:,2,15,:) = SPREAD( zw, 2, statistic_regions+1 )
1743
1744          CASE ( 'w"pt"' )
1745             dopr_index(i) = 16
1746             dopr_unit(i)  = 'K m/s'
1747             hom(:,2,16,:) = SPREAD( zw, 2, statistic_regions+1 )
1748
1749          CASE ( 'w*pt*' )
1750             dopr_index(i) = 17
1751             dopr_unit(i)  = 'K m/s'
1752             hom(:,2,17,:) = SPREAD( zw, 2, statistic_regions+1 )
1753
1754          CASE ( 'wpt' )
1755             dopr_index(i) = 18
1756             dopr_unit(i)  = 'K m/s'
1757             hom(:,2,18,:) = SPREAD( zw, 2, statistic_regions+1 )
1758
1759          CASE ( 'wu' )
1760             dopr_index(i) = 19
1761             dopr_unit(i)  = 'm2/s2'
1762             hom(:,2,19,:) = SPREAD( zw, 2, statistic_regions+1 )
1763             IF ( prandtl_layer )  hom(nzb,2,19,:) = zu(1)
1764
1765          CASE ( 'wv' )
1766             dopr_index(i) = 20
1767             dopr_unit(i)  = 'm2/s2'
1768             hom(:,2,20,:) = SPREAD( zw, 2, statistic_regions+1 )
1769             IF ( prandtl_layer )  hom(nzb,2,20,:) = zu(1)
1770
1771          CASE ( 'w*pt*BC' )
1772             dopr_index(i) = 21
1773             dopr_unit(i)  = 'K m/s'
1774             hom(:,2,21,:) = SPREAD( zw, 2, statistic_regions+1 )
1775
1776          CASE ( 'wptBC' )
1777             dopr_index(i) = 22
1778             dopr_unit(i)  = 'K m/s'
1779             hom(:,2,22,:) = SPREAD( zw, 2, statistic_regions+1 )
1780
1781          CASE ( 'sa', '#sa' )
1782             IF ( .NOT. ocean )  THEN
1783                IF ( myid == 0 )  THEN
1784                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1785                           data_output_pr(i),                          &
1786                           '    is not implemented for ocean = FALSE'
1787                ENDIF
1788                CALL local_stop
1789             ELSE
1790                dopr_index(i) = 23
1791                dopr_unit(i)  = 'psu'
1792                hom(:,2,23,:) = SPREAD( zu, 2, statistic_regions+1 )
1793                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1794                   dopr_initial_index(i) = 26
1795                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1796                   hom(nzb,2,26,:)       = 0.0    ! weil zu(nzb) negativ ist
1797                   data_output_pr(i)     = data_output_pr(i)(2:)
1798                ENDIF
1799             ENDIF
1800
1801          CASE ( 'u*2' )
1802             dopr_index(i) = 30
1803             dopr_unit(i)  = 'm2/s2'
1804             hom(:,2,30,:) = SPREAD( zu, 2, statistic_regions+1 )
1805
1806          CASE ( 'v*2' )
1807             dopr_index(i) = 31
1808             dopr_unit(i)  = 'm2/s2'
1809             hom(:,2,31,:) = SPREAD( zu, 2, statistic_regions+1 )
1810
1811          CASE ( 'w*2' )
1812             dopr_index(i) = 32
1813             dopr_unit(i)  = 'm2/s2'
1814             hom(:,2,32,:) = SPREAD( zw, 2, statistic_regions+1 )
1815
1816          CASE ( 'pt*2' )
1817             dopr_index(i) = 33
1818             dopr_unit(i)  = 'K2'
1819             hom(:,2,33,:) = SPREAD( zu, 2, statistic_regions+1 )
1820
1821          CASE ( 'e*' )
1822             dopr_index(i) = 34
1823             dopr_unit(i)  = 'm2/s2'
1824             hom(:,2,34,:) = SPREAD( zu, 2, statistic_regions+1 )
1825
1826          CASE ( 'w*2pt*' )
1827             dopr_index(i) = 35
1828             dopr_unit(i)  = 'K m2/s2'
1829             hom(:,2,35,:) = SPREAD( zw, 2, statistic_regions+1 )
1830
1831          CASE ( 'w*pt*2' )
1832             dopr_index(i) = 36
1833             dopr_unit(i)  = 'K2 m/s'
1834             hom(:,2,36,:) = SPREAD( zw, 2, statistic_regions+1 )
1835
1836          CASE ( 'w*e*' )
1837             dopr_index(i) = 37
1838             dopr_unit(i)  = 'm3/s3'
1839             hom(:,2,37,:) = SPREAD( zw, 2, statistic_regions+1 )
1840
1841          CASE ( 'w*3' )
1842             dopr_index(i) = 38
1843             dopr_unit(i)  = 'm3/s3'
1844             hom(:,2,38,:) = SPREAD( zw, 2, statistic_regions+1 )
1845
1846          CASE ( 'Sw' )
1847             dopr_index(i) = 39
1848             dopr_unit(i)  = 'none'
1849             hom(:,2,39,:) = SPREAD( zw, 2, statistic_regions+1 )
1850
1851          CASE ( 'q', '#q' )
1852             IF ( .NOT. humidity )  THEN
1853                IF ( myid == 0 )  THEN
1854                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1855                           data_output_pr(i),                          &
1856                           '    is not implemented for humidity = FALSE'
1857                ENDIF
1858                CALL local_stop
1859             ELSE
1860                dopr_index(i) = 41
1861                dopr_unit(i)  = 'kg/kg'
1862                hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 )
1863                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1864                   dopr_initial_index(i) = 26
1865                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1866                   hom(nzb,2,26,:)       = 0.0    ! weil zu(nzb) negativ ist
1867                   data_output_pr(i)     = data_output_pr(i)(2:)
1868                ENDIF
1869             ENDIF
1870
1871          CASE ( 's', '#s' )
1872             IF ( .NOT. passive_scalar )  THEN
1873                IF ( myid == 0 )  THEN
1874                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1875                           data_output_pr(i),                          &
1876                           '    is not implemented for passive_scalar = FALSE'
1877                ENDIF
1878                CALL local_stop
1879             ELSE
1880                dopr_index(i) = 41
1881                dopr_unit(i)  = 'kg/m3'
1882                hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 )
1883                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1884                   dopr_initial_index(i) = 26
1885                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1886                   hom(nzb,2,26,:)       = 0.0    ! weil zu(nzb) negativ ist
1887                   data_output_pr(i)     = data_output_pr(i)(2:)
1888                ENDIF
1889             ENDIF
1890
1891          CASE ( 'qv', '#qv' )
1892             IF ( .NOT. cloud_physics ) THEN
1893                dopr_index(i) = 41
1894                dopr_unit(i)  = 'kg/kg'
1895                hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 )
1896                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1897                   dopr_initial_index(i) = 26
1898                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1899                   hom(nzb,2,26,:)       = 0.0    ! weil zu(nzb) negativ ist
1900                   data_output_pr(i)     = data_output_pr(i)(2:)
1901                ENDIF
1902             ELSE
1903                dopr_index(i) = 42
1904                dopr_unit(i)  = 'kg/kg'
1905                hom(:,2,42,:) = SPREAD( zu, 2, statistic_regions+1 )
1906                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1907                   dopr_initial_index(i) = 27
1908                   hom(:,2,27,:)         = SPREAD( zu, 2, statistic_regions+1 )
1909                   hom(nzb,2,27,:)       = 0.0    ! weil zu(nzb) negativ ist
1910                   data_output_pr(i)     = data_output_pr(i)(2:)
1911                ENDIF
1912             ENDIF
1913
1914          CASE ( 'lpt', '#lpt' )
1915             IF ( .NOT. cloud_physics ) THEN
1916                IF ( myid == 0 )  THEN
1917                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1918                           data_output_pr(i),                          &
1919                           '    is not implemented for cloud_physics = FALSE'
1920                ENDIF
1921                CALL local_stop
1922             ELSE
1923                dopr_index(i) = 4
1924                dopr_unit(i)  = 'K'
1925                hom(:,2,4,:)  = SPREAD( zu, 2, statistic_regions+1 )
1926                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1927                   dopr_initial_index(i) = 7
1928                   hom(:,2,7,:)          = SPREAD( zu, 2, statistic_regions+1 )
1929                   hom(nzb,2,7,:)        = 0.0    ! weil zu(nzb) negativ ist
1930                   data_output_pr(i)     = data_output_pr(i)(2:)
1931                ENDIF
1932             ENDIF
1933
1934          CASE ( 'vpt', '#vpt' )
1935             dopr_index(i) = 44
1936             dopr_unit(i)  = 'K'
1937             hom(:,2,44,:) = SPREAD( zu, 2, statistic_regions+1 )
1938             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1939                dopr_initial_index(i) = 29
1940                hom(:,2,29,:)         = SPREAD( zu, 2, statistic_regions+1 )
1941                hom(nzb,2,29,:)       = 0.0    ! weil zu(nzb) negativ ist
1942                data_output_pr(i)     = data_output_pr(i)(2:)
1943             ENDIF
1944
1945          CASE ( 'w"vpt"' )
1946             dopr_index(i) = 45
1947             dopr_unit(i)  = 'K m/s'
1948             hom(:,2,45,:) = SPREAD( zw, 2, statistic_regions+1 )
1949
1950          CASE ( 'w*vpt*' )
1951             dopr_index(i) = 46
1952             dopr_unit(i)  = 'K m/s'
1953             hom(:,2,46,:) = SPREAD( zw, 2, statistic_regions+1 )
1954
1955          CASE ( 'wvpt' )
1956             dopr_index(i) = 47
1957             dopr_unit(i)  = 'K m/s'
1958             hom(:,2,47,:) = SPREAD( zw, 2, statistic_regions+1 )
1959
1960          CASE ( 'w"q"' )
1961             IF ( .NOT. humidity )  THEN
1962                IF ( myid == 0 )  THEN
1963                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1964                           data_output_pr(i),                          &
1965                           '    is not implemented for humidity = FALSE'
1966                ENDIF
1967                CALL local_stop
1968             ELSE
1969                dopr_index(i) = 48
1970                dopr_unit(i)  = 'kg/kg m/s'
1971                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
1972             ENDIF
1973
1974          CASE ( 'w*q*' )
1975             IF ( .NOT. humidity )  THEN
1976                IF ( myid == 0 )  THEN
1977                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1978                           data_output_pr(i),                          &
1979                           '    is not implemented for humidity = FALSE'
1980                ENDIF
1981                CALL local_stop
1982             ELSE
1983                dopr_index(i) = 49
1984                dopr_unit(i)  = 'kg/kg m/s'
1985                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
1986             ENDIF
1987
1988          CASE ( 'wq' )
1989             IF ( .NOT. humidity )  THEN
1990                IF ( myid == 0 )  THEN
1991                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1992                           data_output_pr(i),                          &
1993                           '    is not implemented for humidity = FALSE'
1994                ENDIF
1995                CALL local_stop
1996             ELSE
1997                dopr_index(i) = 50
1998                dopr_unit(i)  = 'kg/kg m/s'
1999                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
2000             ENDIF
2001
2002          CASE ( 'w"s"' )
2003             IF ( .NOT. passive_scalar ) THEN
2004                IF ( myid == 0 )  THEN
2005                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2006                           data_output_pr(i),                          &
2007                           '    is not implemented for passive_scalar = FALSE'
2008                ENDIF
2009                CALL local_stop
2010             ELSE
2011                dopr_index(i) = 48
2012                dopr_unit(i)  = 'kg/m3 m/s'
2013                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
2014             ENDIF
2015
2016          CASE ( 'w*s*' )
2017             IF ( .NOT. passive_scalar ) THEN
2018                IF ( myid == 0 )  THEN
2019                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2020                           data_output_pr(i),                          &
2021                           '    is not implemented for passive_scalar = FALSE'
2022                ENDIF
2023                CALL local_stop
2024             ELSE
2025                dopr_index(i) = 49
2026                dopr_unit(i)  = 'kg/m3 m/s'
2027                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
2028             ENDIF
2029
2030          CASE ( 'ws' )
2031             IF ( .NOT. passive_scalar ) THEN
2032                IF ( myid == 0 )  THEN
2033                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2034                           data_output_pr(i),                          &
2035                           '    is not implemented for passive_scalar = FALSE'
2036                ENDIF
2037                CALL local_stop
2038             ELSE
2039                dopr_index(i) = 50
2040                dopr_unit(i)  = 'kg/m3 m/s'
2041                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
2042             ENDIF
2043
2044          CASE ( 'w"qv"' )
2045             IF ( humidity  .AND.  .NOT. cloud_physics ) &
2046             THEN
2047                dopr_index(i) = 48
2048                dopr_unit(i)  = 'kg/kg m/s'
2049                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
2050             ELSEIF( humidity .AND. cloud_physics ) THEN
2051                dopr_index(i) = 51
2052                dopr_unit(i)  = 'kg/kg m/s'
2053                hom(:,2,51,:) = SPREAD( zw, 2, statistic_regions+1 )
2054             ELSE
2055                IF ( myid == 0 )  THEN
2056                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2057                           data_output_pr(i),                          &
2058                           '    is not implemented for cloud_physics = FALSE', &
2059                           '    and                    humidity      = FALSE'
2060                ENDIF
2061                CALL local_stop                   
2062             ENDIF
2063
2064          CASE ( 'w*qv*' )
2065             IF ( humidity  .AND.  .NOT. cloud_physics ) &
2066             THEN
2067                dopr_index(i) = 49
2068                dopr_unit(i)  = 'kg/kg m/s'
2069                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
2070             ELSEIF( humidity .AND. cloud_physics ) THEN
2071                dopr_index(i) = 52
2072                dopr_unit(i)  = 'kg/kg m/s'
2073                hom(:,2,52,:) = SPREAD( zw, 2, statistic_regions+1 )
2074             ELSE
2075                IF ( myid == 0 )  THEN
2076                   PRINT*, '+++ check_parameters:  data_output_pr = ',         &
2077                           data_output_pr(i),                                  &
2078                           '    is not implemented for cloud_physics = FALSE', &
2079                           '                       and humidity      = FALSE'
2080                ENDIF
2081                CALL local_stop                   
2082             ENDIF
2083
2084          CASE ( 'wqv' )
2085             IF ( humidity  .AND.  .NOT. cloud_physics ) &
2086             THEN
2087                dopr_index(i) = 50
2088                dopr_unit(i)  = 'kg/kg m/s'
2089                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
2090             ELSEIF( humidity .AND. cloud_physics ) THEN
2091                dopr_index(i) = 53
2092                dopr_unit(i)  = 'kg/kg m/s'
2093                hom(:,2,53,:) = SPREAD( zw, 2, statistic_regions+1 )
2094             ELSE
2095                IF ( myid == 0 )  THEN
2096                   PRINT*, '+++ check_parameters:  data_output_pr = ',         &
2097                           data_output_pr(i),                                  &
2098                           '    is not implemented for cloud_physics = FALSE', &
2099                           '                       and humidity      = FALSE'
2100                ENDIF
2101                CALL local_stop                   
2102             ENDIF
2103
2104          CASE ( 'ql' )
2105             IF ( .NOT. cloud_physics  .AND.  .NOT. cloud_droplets )  THEN
2106                IF ( myid == 0 )  THEN
2107                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2108                           data_output_pr(i),                          &
2109                           '    is not implemented for cloud_physics = FALSE'
2110                ENDIF
2111                CALL local_stop
2112             ELSE
2113                dopr_index(i) = 54
2114                dopr_unit(i)  = 'kg/kg'
2115                hom(:,2,54,:)  = SPREAD( zu, 2, statistic_regions+1 )
2116             ENDIF
2117
2118          CASE ( 'w*u*u*/dz' )
2119             dopr_index(i) = 55
2120             dopr_unit(i)  = 'm2/s3'
2121             hom(:,2,55,:) = SPREAD( zu, 2, statistic_regions+1 )
2122
2123          CASE ( 'w*p*/dz' )
2124             dopr_index(i) = 56
2125             dopr_unit(i)  = 'm2/s3'
2126             hom(:,2,56,:) = SPREAD( zw, 2, statistic_regions+1 )
2127
2128          CASE ( 'w"e/dz' )
2129             dopr_index(i) = 57
2130             dopr_unit(i)  = 'm2/s3'
2131             hom(:,2,57,:) = SPREAD( zu, 2, statistic_regions+1 )
2132
2133          CASE ( 'u"pt"' )
2134             dopr_index(i) = 58
2135             dopr_unit(i)  = 'K m/s'
2136             hom(:,2,58,:) = SPREAD( zu, 2, statistic_regions+1 )
2137
2138          CASE ( 'u*pt*' )
2139             dopr_index(i) = 59
2140             dopr_unit(i)  = 'K m/s'
2141             hom(:,2,59,:) = SPREAD( zu, 2, statistic_regions+1 )
2142
2143          CASE ( 'upt_t' )
2144             dopr_index(i) = 60
2145             dopr_unit(i)  = 'K m/s'
2146             hom(:,2,60,:) = SPREAD( zu, 2, statistic_regions+1 )
2147
2148          CASE ( 'v"pt"' )
2149             dopr_index(i) = 61
2150             dopr_unit(i)  = 'K m/s'
2151             hom(:,2,61,:) = SPREAD( zu, 2, statistic_regions+1 )
2152             
2153          CASE ( 'v*pt*' )
2154             dopr_index(i) = 62
2155             dopr_unit(i)  = 'K m/s'
2156             hom(:,2,62,:) = SPREAD( zu, 2, statistic_regions+1 )
2157
2158          CASE ( 'vpt_t' )
2159             dopr_index(i) = 63
2160             dopr_unit(i)  = 'K m/s'
2161             hom(:,2,63,:) = SPREAD( zu, 2, statistic_regions+1 )
2162
2163          CASE ( 'rho' )
2164             dopr_index(i) = 64
2165             dopr_unit(i)  = 'kg/m3'
2166             hom(:,2,64,:) = SPREAD( zu, 2, statistic_regions+1 )
2167
2168          CASE ( 'w"sa"' )
2169             IF ( .NOT. ocean ) THEN
2170                IF ( myid == 0 )  THEN
2171                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2172                           data_output_pr(i),                          &
2173                           '    is not implemented for ocean = FALSE'
2174                ENDIF
2175                CALL local_stop
2176             ELSE
2177                dopr_index(i) = 65
2178                dopr_unit(i)  = 'psu m/s'
2179                hom(:,2,65,:) = SPREAD( zw, 2, statistic_regions+1 )
2180             ENDIF
2181
2182          CASE ( 'w*sa*' )
2183             IF ( .NOT. ocean ) THEN
2184                IF ( myid == 0 )  THEN
2185                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2186                           data_output_pr(i),                          &
2187                           '    is not implemented for ocean = FALSE'
2188                ENDIF
2189                CALL local_stop
2190             ELSE
2191                dopr_index(i) = 66
2192                dopr_unit(i)  = 'psu m/s'
2193                hom(:,2,66,:) = SPREAD( zw, 2, statistic_regions+1 )
2194             ENDIF
2195
2196          CASE ( 'wsa' )
2197             IF ( .NOT. ocean ) THEN
2198                IF ( myid == 0 )  THEN
2199                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
2200                           data_output_pr(i),                          &
2201                           '    is not implemented for ocean = FALSE'
2202                ENDIF
2203                CALL local_stop
2204             ELSE
2205                dopr_index(i) = 67
2206                dopr_unit(i)  = 'psu m/s'
2207                hom(:,2,67,:) = SPREAD( zw, 2, statistic_regions+1 )
2208             ENDIF
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
2221          CASE DEFAULT
2222
2223             CALL user_check_data_output_pr( data_output_pr(i), i, unit )
2224
2225             IF ( unit == 'illegal' )  THEN
2226                IF ( myid == 0 )  THEN
2227                   IF ( data_output_pr_user(1) /= ' ' )  THEN
2228                      PRINT*, '+++ check_parameters:  illegal value for data_',&
2229                                   'output_pr or data_output_pr_user: "',      &
2230                                   TRIM( data_output_pr(i) ), '"'
2231                   ELSE
2232                      PRINT*, '+++ check_parameters:  illegal value for data_',&
2233                                   'output_pr: "', TRIM( data_output_pr(i) ),'"'
2234                   ENDIF
2235                ENDIF
2236                CALL local_stop
2237             ENDIF
2238
2239       END SELECT
2240!
2241!--    Check to which of the predefined coordinate systems the profile belongs
2242       DO  k = 1, crmax
2243          IF ( INDEX( cross_profiles(k), ' '//TRIM( data_output_pr(i) )//' ' ) &
2244               /=0 ) &
2245          THEN
2246             dopr_crossindex(i) = k
2247             EXIT
2248          ENDIF
2249       ENDDO
2250!
2251!--    Generate the text for the labels of the PROFIL output file. "-characters
2252!--    must be substituted, otherwise PROFIL would interpret them as TeX
2253!--    control characters
2254       dopr_label(i) = data_output_pr(i)
2255       position = INDEX( dopr_label(i) , '"' )
2256       DO WHILE ( position /= 0 )
2257          dopr_label(i)(position:position) = ''''
2258          position = INDEX( dopr_label(i) , '"' )
2259       ENDDO
2260
2261    ENDDO
2262
2263!
2264!-- y-value range of the coordinate system (PROFIL).
2265!-- x-value range determined in plot_1d.
2266    IF ( .NOT. ocean )  THEN
2267       cross_uymin = 0.0
2268       IF ( z_max_do1d == -1.0 )  THEN
2269          cross_uymax = zu(nzt+1)
2270       ELSEIF ( z_max_do1d < zu(nzb+1)  .OR.  z_max_do1d > zu(nzt+1) )  THEN
2271          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  z_max_do1d=',  &
2272                                    z_max_do1d,' must be >= ', zu(nzb+1),  &
2273                                    ' or <= ', zu(nzt+1)
2274          CALL local_stop
2275       ELSE
2276          cross_uymax = z_max_do1d
2277       ENDIF
2278    ENDIF
2279
2280!
2281!-- Check whether the chosen normalizing factor for the coordinate systems is
2282!-- permissible
2283    DO  i = 1, crmax
2284       SELECT CASE ( TRIM( cross_normalized_x(i) ) )  ! TRIM required on IBM
2285
2286          CASE ( '', 'wpt0', 'ws2', 'tsw2', 'ws3', 'ws2tsw', 'wstsw2' )
2287             j = 0
2288
2289          CASE DEFAULT
2290             IF ( myid == 0 )  THEN
2291                PRINT*, '+++ check_parameters: unknown normalize method'
2292                PRINT*, '    cross_normalized_x="',cross_normalized_x(i),'"'
2293             ENDIF
2294             CALL local_stop
2295
2296       END SELECT
2297       SELECT CASE ( TRIM( cross_normalized_y(i) ) )  ! TRIM required on IBM
2298
2299          CASE ( '', 'z_i' )
2300             j = 0
2301
2302          CASE DEFAULT
2303             IF ( myid == 0 )  THEN
2304                PRINT*, '+++ check_parameters: unknown normalize method'
2305                PRINT*, '    cross_normalized_y="',cross_normalized_y(i),'"'
2306             ENDIF
2307             CALL local_stop
2308
2309       END SELECT
2310    ENDDO
2311!
2312!-- Check normalized y-value range of the coordinate system (PROFIL)
2313    IF ( z_max_do1d_normalized /= -1.0  .AND.  z_max_do1d_normalized <= 0.0 ) &
2314    THEN
2315       IF ( myid == 0 )  PRINT*,'+++ check_parameters:  z_max_do1d_normalize', &
2316                                'd=', z_max_do1d_normalized, ' must be >= 0.0'
2317       CALL local_stop
2318    ENDIF
2319
2320
2321!
2322!-- Append user-defined data output variables to the standard data output
2323    IF ( data_output_user(1) /= ' ' )  THEN
2324       i = 1
2325       DO  WHILE ( data_output(i) /= ' '  .AND.  i <= 100 )
2326          i = i + 1
2327       ENDDO
2328       j = 1
2329       DO  WHILE ( data_output_user(j) /= ' '  .AND.  j <= 100 )
2330          IF ( i > 100 )  THEN
2331             IF ( myid == 0 )  THEN
2332                PRINT*, '+++ check_parameters: number of output quantitities', &
2333                             ' given by data_output and data_output_user'
2334                PRINT*, '                      exceeds the limit of 100'
2335             ENDIF
2336             CALL local_stop
2337          ENDIF
2338          data_output(i) = data_output_user(j)
2339          i = i + 1
2340          j = j + 1
2341       ENDDO
2342    ENDIF
2343
2344!
2345!-- Check and set steering parameters for 2d/3d data output and averaging
2346    i   = 1
2347    DO  WHILE ( data_output(i) /= ' '  .AND.  i <= 100 )
2348!
2349!--    Check for data averaging
2350       ilen = LEN_TRIM( data_output(i) )
2351       j = 0                                                 ! no data averaging
2352       IF ( ilen > 3 )  THEN
2353          IF ( data_output(i)(ilen-2:ilen) == '_av' )  THEN
2354             j = 1                                           ! data averaging
2355             data_output(i) = data_output(i)(1:ilen-3)
2356          ENDIF
2357       ENDIF
2358!
2359!--    Check for cross section or volume data
2360       ilen = LEN_TRIM( data_output(i) )
2361       k = 0                                                   ! 3d data
2362       var = data_output(i)(1:ilen)
2363       IF ( ilen > 3 )  THEN
2364          IF ( data_output(i)(ilen-2:ilen) == '_xy'  .OR. &
2365               data_output(i)(ilen-2:ilen) == '_xz'  .OR. &
2366               data_output(i)(ilen-2:ilen) == '_yz' )  THEN
2367             k = 1                                             ! 2d data
2368             var = data_output(i)(1:ilen-3)
2369          ENDIF
2370       ENDIF
2371!
2372!--    Check for allowed value and set units
2373       SELECT CASE ( TRIM( var ) )
2374
2375          CASE ( 'e' )
2376             IF ( constant_diffusion )  THEN
2377                IF ( myid == 0 )  THEN
2378                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2379                                '" requires constant_diffusion = .FALSE.'
2380                ENDIF
2381                CALL local_stop
2382             ENDIF
2383             unit = 'm2/s2'
2384
2385          CASE ( 'pc', 'pr' )
2386             IF ( .NOT. particle_advection )  THEN
2387                IF ( myid == 0 )  THEN
2388                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2389                                '" requires a "particles_par"-NAMELIST'
2390                   PRINT*, '                      in the parameter file (PARIN)'
2391                ENDIF
2392                CALL local_stop
2393             ENDIF
2394             IF ( TRIM( var ) == 'pc' )  unit = 'number'
2395             IF ( TRIM( var ) == 'pr' )  unit = 'm'
2396
2397          CASE ( 'q', 'vpt' )
2398             IF ( .NOT. humidity )  THEN
2399                IF ( myid == 0 )  THEN
2400                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2401                                '" requires humidity = .TRUE.'
2402                ENDIF
2403                CALL local_stop
2404             ENDIF
2405             IF ( TRIM( var ) == 'q'   )  unit = 'kg/kg'
2406             IF ( TRIM( var ) == 'vpt' )  unit = 'K'
2407
2408          CASE ( 'ql' )
2409             IF ( .NOT. ( cloud_physics  .OR.  cloud_droplets ) )  THEN
2410                IF ( myid == 0 )  THEN
2411                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2412                                '" requires cloud_physics = .TRUE.'
2413                   PRINT*, '                      or cloud_droplets = .TRUE.'
2414                ENDIF
2415                CALL local_stop
2416             ENDIF
2417             unit = 'kg/kg'
2418
2419          CASE ( 'ql_c', 'ql_v', 'ql_vp' )
2420             IF ( .NOT. cloud_droplets )  THEN
2421                IF ( myid == 0 )  THEN
2422                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2423                                '" requires cloud_droplets = .TRUE.'
2424                ENDIF
2425                CALL local_stop
2426             ENDIF
2427             IF ( TRIM( var ) == 'ql_c'  )  unit = 'kg/kg'
2428             IF ( TRIM( var ) == 'ql_v'  )  unit = 'm3'
2429             IF ( TRIM( var ) == 'ql_vp' )  unit = 'none'
2430
2431          CASE ( 'qv' )
2432             IF ( .NOT. cloud_physics )  THEN
2433                IF ( myid == 0 )  THEN
2434                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2435                                '" requires cloud_physics = .TRUE.'
2436                ENDIF
2437                CALL local_stop
2438             ENDIF
2439             unit = 'kg/kg'
2440
2441          CASE ( 'rho' )
2442             IF ( .NOT. ocean )  THEN
2443                IF ( myid == 0 )  THEN
2444                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2445                                '" requires ocean = .TRUE.'
2446                ENDIF
2447                CALL local_stop
2448             ENDIF
2449             unit = 'kg/m3'
2450
2451          CASE ( 's' )
2452             IF ( .NOT. passive_scalar )  THEN
2453                IF ( myid == 0 )  THEN
2454                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2455                                '" requires passive_scalar = .TRUE.'
2456                ENDIF
2457                CALL local_stop
2458             ENDIF
2459             unit = 'conc'
2460
2461          CASE ( 'sa' )
2462             IF ( .NOT. ocean )  THEN
2463                IF ( myid == 0 )  THEN
2464                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2465                                '" requires ocean = .TRUE.'
2466                ENDIF
2467                CALL local_stop
2468             ENDIF
2469             unit = 'psu'
2470
2471          CASE ( 'u*', 't*', 'lwp*', 'pra*', 'prr*', 'z0*' )
2472             IF ( k == 0  .OR.  data_output(i)(ilen-2:ilen) /= '_xy' )  THEN
2473                IF ( myid == 0 )  THEN
2474                   PRINT*, '+++ check_parameters:  illegal value for data_',&
2475                                'output: "', TRIM( var ), '" is only allowed'
2476                   PRINT*, '                       for horizontal cross section'
2477                ENDIF
2478                CALL local_stop
2479             ENDIF
2480             IF ( TRIM( var ) == 'lwp*'  .AND.  .NOT. cloud_physics )  THEN
2481                IF ( myid == 0 )  THEN
2482                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2483                                '" requires cloud_physics = .TRUE.'
2484                ENDIF
2485                CALL local_stop
2486             ENDIF
2487             IF ( TRIM( var ) == 'pra*'  .AND.  .NOT. precipitation )  THEN
2488                IF ( myid == 0 )  THEN
2489                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2490                                '" requires precipitation = .TRUE.'
2491                ENDIF
2492                CALL local_stop
2493             ENDIF
2494             IF ( TRIM( var ) == 'pra*'  .AND.  j == 1 )  THEN
2495                IF ( myid == 0 )  THEN
2496                   PRINT*, '+++ check_parameters: temporal averaging of ', &
2497                           ' precipitation amount "', TRIM( var ),         &
2498                           '" not possible' 
2499                ENDIF
2500                CALL local_stop
2501             ENDIF
2502             IF ( TRIM( var ) == 'prr*'  .AND.  .NOT. precipitation )  THEN
2503                IF ( myid == 0 )  THEN
2504                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
2505                                '" requires precipitation = .TRUE.'
2506                ENDIF
2507                CALL local_stop
2508             ENDIF
2509
2510
2511             IF ( TRIM( var ) == 'u*'   )  unit = 'm/s'
2512             IF ( TRIM( var ) == 't*'   )  unit = 'K'
2513             IF ( TRIM( var ) == 'lwp*' )  unit = 'kg/kg*m'
2514             IF ( TRIM( var ) == 'pra*' )  unit = 'mm'
2515             IF ( TRIM( var ) == 'prr*' )  unit = 'mm/s'
2516             IF ( TRIM( var ) == 'z0*'  )  unit = 'm'
2517
2518          CASE ( 'p', 'pt', 'u', 'v', 'w' )
2519             IF ( TRIM( var ) == 'p'  )  unit = 'Pa'
2520             IF ( TRIM( var ) == 'pt' )  unit = 'K'
2521             IF ( TRIM( var ) == 'u'  )  unit = 'm/s'
2522             IF ( TRIM( var ) == 'v'  )  unit = 'm/s'
2523             IF ( TRIM( var ) == 'w'  )  unit = 'm/s'
2524             CONTINUE
2525
2526          CASE DEFAULT
2527             CALL user_check_data_output( var, unit )
2528
2529             IF ( unit == 'illegal' )  THEN
2530                IF ( myid == 0 )  THEN
2531                   IF ( data_output_user(1) /= ' ' )  THEN
2532                      PRINT*, '+++ check_parameters:  illegal value for data_',&
2533                                   'output or data_output_user: "',            &
2534                                   TRIM( data_output(i) ), '"'
2535                   ELSE
2536                      PRINT*, '+++ check_parameters:  illegal value for data_',&
2537                                   'output: "', TRIM( data_output(i) ), '"'
2538                   ENDIF
2539                ENDIF
2540                CALL local_stop
2541             ENDIF
2542
2543       END SELECT
2544!
2545!--    Set the internal steering parameters appropriately
2546       IF ( k == 0 )  THEN
2547          do3d_no(j)              = do3d_no(j) + 1
2548          do3d(j,do3d_no(j))      = data_output(i)
2549          do3d_unit(j,do3d_no(j)) = unit
2550       ELSE
2551          do2d_no(j)              = do2d_no(j) + 1
2552          do2d(j,do2d_no(j))      = data_output(i)
2553          do2d_unit(j,do2d_no(j)) = unit
2554          IF ( data_output(i)(ilen-2:ilen) == '_xy' )  THEN
2555             data_output_xy(j) = .TRUE.
2556          ENDIF
2557          IF ( data_output(i)(ilen-2:ilen) == '_xz' )  THEN
2558             data_output_xz(j) = .TRUE.
2559          ENDIF
2560          IF ( data_output(i)(ilen-2:ilen) == '_yz' )  THEN
2561             data_output_yz(j) = .TRUE.
2562          ENDIF
2563       ENDIF
2564
2565       IF ( j == 1 )  THEN
2566!
2567!--       Check, if variable is already subject to averaging
2568          found = .FALSE.
2569          DO  k = 1, doav_n
2570             IF ( TRIM( doav(k) ) == TRIM( var ) )  found = .TRUE.
2571          ENDDO
2572
2573          IF ( .NOT. found )  THEN
2574             doav_n = doav_n + 1
2575             doav(doav_n) = var
2576          ENDIF
2577       ENDIF
2578
2579       i = i + 1
2580    ENDDO
2581
2582!
2583!-- Store sectional planes in one shared array
2584    section(:,1) = section_xy
2585    section(:,2) = section_xz
2586    section(:,3) = section_yz
2587
2588!
2589!-- Upper plot limit (grid point value) for 1D profiles
2590    IF ( z_max_do1d == -1.0 )  THEN
2591       nz_do1d = nzt+1
2592    ELSE
2593       DO  k = nzb+1, nzt+1
2594          nz_do1d = k
2595          IF ( zw(k) > z_max_do1d )  EXIT
2596       ENDDO
2597    ENDIF
2598
2599!
2600!-- Upper plot limit for 2D vertical sections
2601    IF ( z_max_do2d == -1.0 )  z_max_do2d = zu(nzt)
2602    IF ( z_max_do2d < zu(nzb+1)  .OR.  z_max_do2d > zu(nzt) )  THEN
2603       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  z_max_do2d=',        &
2604                                 z_max_do2d, ' must be >= ', zu(nzb+1),       &
2605                                 '(zu(nzb+1)) and <= ', zu(nzt), ' (zu(nzt))'
2606       CALL local_stop
2607    ENDIF
2608
2609!
2610!-- Upper plot limit for 3D arrays
2611    IF ( nz_do3d == -9999 )  nz_do3d = nzt + 1
2612
2613!
2614!-- Determine and check accuracy for compressed 3D plot output
2615    IF ( do3d_compress )  THEN
2616!
2617!--    Compression only permissible on T3E machines
2618       IF ( host(1:3) /= 't3e' )  THEN
2619          IF ( myid == 0 )  THEN
2620             PRINT*, '+++ check_parameters: do3d_compress = .TRUE. not allow', &
2621                          'ed on host "', TRIM( host ), '"'
2622          ENDIF
2623          CALL local_stop
2624       ENDIF
2625
2626       i = 1
2627       DO  WHILE ( do3d_comp_prec(i) /= ' ' )
2628
2629          ilen = LEN_TRIM( do3d_comp_prec(i) )
2630          IF ( LLT( do3d_comp_prec(i)(ilen:ilen), '0' ) .OR. &
2631               LGT( do3d_comp_prec(i)(ilen:ilen), '9' ) )  THEN
2632             IF ( myid == 0 )  THEN
2633                PRINT*, '+++ check_parameters: illegal precision: ', &
2634                        'do3d_comp_prec(', i, ')="', TRIM(do3d_comp_prec(i)),'"'
2635             ENDIF
2636             CALL local_stop
2637          ENDIF
2638
2639          prec = IACHAR( do3d_comp_prec(i)(ilen:ilen) ) - IACHAR( '0' )
2640          var = do3d_comp_prec(i)(1:ilen-1)
2641
2642          SELECT CASE ( var )
2643
2644             CASE ( 'u' )
2645                j = 1
2646             CASE ( 'v' )
2647                j = 2
2648             CASE ( 'w' )
2649                j = 3
2650             CASE ( 'p' )
2651                j = 4
2652             CASE ( 'pt' )
2653                j = 5
2654
2655             CASE DEFAULT
2656                IF ( myid == 0 )  THEN
2657                   PRINT*, '+++ check_parameters: unknown variable in ', &
2658                               'assignment'
2659                   PRINT*, '    do3d_comp_prec(', i, ')="', &
2660                           TRIM( do3d_comp_prec(i) ),'"'
2661                ENDIF
2662                CALL local_stop               
2663
2664          END SELECT
2665
2666          plot_3d_precision(j)%precision = prec
2667          i = i + 1
2668
2669       ENDDO
2670    ENDIF
2671
2672!
2673!-- Check the data output format(s)
2674    IF ( data_output_format(1) == ' ' )  THEN
2675!
2676!--    Default value
2677       netcdf_output = .TRUE.
2678    ELSE
2679       i = 1
2680       DO  WHILE ( data_output_format(i) /= ' ' )
2681
2682          SELECT CASE ( data_output_format(i) )
2683
2684             CASE ( 'netcdf' )
2685                netcdf_output = .TRUE.
2686             CASE ( 'iso2d' )
2687                iso2d_output  = .TRUE.
2688             CASE ( 'profil' )
2689                profil_output = .TRUE.
2690             CASE ( 'avs' )
2691                avs_output    = .TRUE.
2692
2693             CASE DEFAULT
2694                IF ( myid == 0 )  THEN
2695                   PRINT*, '+++ check_parameters:'
2696                   PRINT*, '    unknown value for data_output_format "', &
2697                                TRIM( data_output_format(i) ),'"'
2698                ENDIF
2699                CALL local_stop               
2700
2701          END SELECT
2702
2703          i = i + 1
2704          IF ( i > 10 )  EXIT
2705
2706       ENDDO
2707
2708    ENDIF
2709
2710!
2711!-- Check netcdf precison
2712    ldum = .FALSE.
2713    CALL define_netcdf_header( 'ch', ldum, 0 )
2714
2715!
2716!-- Check, whether a constant diffusion coefficient shall be used
2717    IF ( km_constant /= -1.0 )  THEN
2718       IF ( km_constant < 0.0 )  THEN
2719          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  km_constant=', &
2720                                    km_constant, ' < 0.0'
2721          CALL local_stop
2722       ELSE
2723          IF ( prandtl_number < 0.0 )  THEN
2724             IF ( myid == 0 )  PRINT*,'+++ check_parameters:  prandtl_number=',&
2725                                      prandtl_number, ' < 0.0'
2726             CALL local_stop
2727          ENDIF
2728          constant_diffusion = .TRUE.
2729
2730          IF ( prandtl_layer )  THEN
2731             IF ( myid == 0 )  PRINT*, '+++ check_parameters:  prandtl_layer ',&
2732                                       'is not allowed with fixed value of km'
2733             CALL local_stop
2734          ENDIF
2735       ENDIF
2736    ENDIF
2737
2738!
2739!-- In case of non-cyclic lateral boundaries, set the default maximum value
2740!-- for the horizontal diffusivity used within the outflow damping layer,
2741!-- and check/set the width of the damping layer
2742    IF ( bc_lr /= 'cyclic' ) THEN
2743       IF ( km_damp_max == -1.0 )  THEN
2744          km_damp_max = 0.5 * dx
2745       ENDIF
2746       IF ( outflow_damping_width == -1.0 )  THEN
2747          outflow_damping_width = MIN( 20, nx/2 )
2748       ENDIF
2749       IF ( outflow_damping_width <= 0  .OR.  outflow_damping_width > nx )  THEN
2750          IF ( myid == 0 )  PRINT*, '+++ check_parameters: outflow_damping w',&
2751                                    'idth out of range'
2752          CALL local_stop
2753       ENDIF
2754    ENDIF
2755
2756    IF ( bc_ns /= 'cyclic' )  THEN
2757       IF ( km_damp_max == -1.0 )  THEN
2758          km_damp_max = 0.5 * dy
2759       ENDIF
2760       IF ( outflow_damping_width == -1.0 )  THEN
2761          outflow_damping_width = MIN( 20, ny/2 )
2762       ENDIF
2763       IF ( outflow_damping_width <= 0  .OR.  outflow_damping_width > ny )  THEN
2764          IF ( myid == 0 )  PRINT*, '+++ check_parameters: outflow_damping w',&
2765                                    'idth out of range'
2766          CALL local_stop
2767       ENDIF
2768    ENDIF
2769
2770!
2771!-- Check value range for rif
2772    IF ( rif_min >= rif_max )  THEN
2773       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  rif_min=', rif_min, &
2774                                 ' must be less than rif_max=', rif_max
2775       CALL local_stop
2776    ENDIF
2777
2778!
2779!-- Determine upper and lower hight level indices for random perturbations
2780    IF ( disturbance_level_b == -9999999.9 )  THEN
2781       IF ( ocean ) THEN
2782          disturbance_level_b     = zu((nzt*2)/3)
2783          disturbance_level_ind_b = ( nzt * 2 ) / 3
2784       ELSE
2785          disturbance_level_b     = zu(nzb+3)
2786          disturbance_level_ind_b = nzb + 3
2787       ENDIF
2788    ELSEIF ( disturbance_level_b < zu(3) )  THEN
2789       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_b=',&
2790                                 disturbance_level_b, ' must be >= ',zu(3),    &
2791                                 '(zu(3))'
2792       CALL local_stop
2793    ELSEIF ( disturbance_level_b > zu(nzt-2) )  THEN
2794       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_b=',&
2795                                 disturbance_level_b, ' must be <= ',zu(nzt-2),&
2796                                 '(zu(nzt-2))'
2797       CALL local_stop
2798    ELSE
2799       DO  k = 3, nzt-2
2800          IF ( disturbance_level_b <= zu(k) )  THEN
2801             disturbance_level_ind_b = k
2802             EXIT
2803          ENDIF
2804       ENDDO
2805    ENDIF
2806
2807    IF ( disturbance_level_t == -9999999.9 )  THEN
2808       IF ( ocean )  THEN
2809          disturbance_level_t     = zu(nzt-3)
2810          disturbance_level_ind_t = nzt - 3
2811       ELSE
2812          disturbance_level_t     = zu(nzt/3)
2813          disturbance_level_ind_t = nzt / 3
2814       ENDIF
2815    ELSEIF ( disturbance_level_t > zu(nzt-2) )  THEN
2816       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_t=',&
2817                                 disturbance_level_t, ' must be <= ',zu(nzt-2),&
2818                                 '(zu(nzt-2))'
2819       CALL local_stop
2820    ELSEIF ( disturbance_level_t < disturbance_level_b )  THEN
2821       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_t=',&
2822                                 disturbance_level_t, ' must be >= ',          &
2823                                 'disturbance_level_b=', disturbance_level_b
2824       CALL local_stop
2825    ELSE
2826       DO  k = 3, nzt-2
2827          IF ( disturbance_level_t <= zu(k) )  THEN
2828             disturbance_level_ind_t = k
2829             EXIT
2830          ENDIF
2831       ENDDO
2832    ENDIF
2833
2834!
2835!-- Check again whether the levels determined this way are ok.
2836!-- Error may occur at automatic determination and too few grid points in
2837!-- z-direction.
2838    IF ( disturbance_level_ind_t < disturbance_level_ind_b )  THEN
2839       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  ',                &
2840                                 'disturbance_level_ind_t=',               &
2841                                 disturbance_level_ind_t, ' must be >= ',  &
2842                                 'disturbance_level_ind_b=',               &
2843                                 disturbance_level_ind_b
2844       CALL local_stop
2845    ENDIF
2846
2847!
2848!-- Determine the horizontal index range for random perturbations.
2849!-- In case of non-cyclic horizontal boundaries, no perturbations are imposed
2850!-- near the inflow and the perturbation area is further limited to ...(1)
2851!-- after the initial phase of the flow.
2852    dist_nxl = 0;  dist_nxr = nx
2853    dist_nys = 0;  dist_nyn = ny
2854    IF ( bc_lr /= 'cyclic' )  THEN
2855       IF ( inflow_disturbance_begin == -1 )  THEN
2856          inflow_disturbance_begin = MIN( 10, nx/2 )
2857       ENDIF
2858       IF ( inflow_disturbance_begin < 0  .OR.  inflow_disturbance_begin > nx )&
2859       THEN
2860          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2861                                    '_begin out of range'
2862          CALL local_stop
2863       ENDIF
2864       IF ( inflow_disturbance_end == -1 )  THEN
2865          inflow_disturbance_end = MIN( 100, 3*nx/4 )
2866       ENDIF
2867       IF ( inflow_disturbance_end < 0  .OR.  inflow_disturbance_end > nx )    &
2868       THEN
2869          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2870                                    '_end out of range'
2871          CALL local_stop
2872       ENDIF
2873    ELSEIF ( bc_ns /= 'cyclic' )  THEN
2874       IF ( inflow_disturbance_begin == -1 )  THEN
2875          inflow_disturbance_begin = MIN( 10, ny/2 )
2876       ENDIF
2877       IF ( inflow_disturbance_begin < 0  .OR.  inflow_disturbance_begin > ny )&
2878       THEN
2879          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2880                                    '_begin out of range'
2881          CALL local_stop
2882       ENDIF
2883       IF ( inflow_disturbance_end == -1 )  THEN
2884          inflow_disturbance_end = MIN( 100, 3*ny/4 )
2885       ENDIF
2886       IF ( inflow_disturbance_end < 0  .OR.  inflow_disturbance_end > ny )    &
2887       THEN
2888          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2889                                    '_end out of range'
2890          CALL local_stop
2891       ENDIF
2892    ENDIF
2893
2894    IF ( bc_lr == 'radiation/dirichlet' )  THEN
2895       dist_nxr    = nx - inflow_disturbance_begin
2896       dist_nxl(1) = nx - inflow_disturbance_end
2897    ELSEIF ( bc_lr == 'dirichlet/radiation' )  THEN
2898       dist_nxl    = inflow_disturbance_begin
2899       dist_nxr(1) = inflow_disturbance_end
2900    ENDIF
2901    IF ( bc_ns == 'dirichlet/radiation' )  THEN
2902       dist_nyn    = ny - inflow_disturbance_begin
2903       dist_nys(1) = ny - inflow_disturbance_end
2904    ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
2905       dist_nys    = inflow_disturbance_begin
2906       dist_nyn(1) = inflow_disturbance_end
2907    ENDIF
2908
2909!
2910!-- Check random generator
2911    IF ( random_generator /= 'system-specific'  .AND. &
2912         random_generator /= 'numerical-recipes' )  THEN
2913       IF ( myid == 0 )  THEN
2914          PRINT*, '+++ check_parameters:'
2915          PRINT*, '    unknown random generator: random_generator=', &
2916                  random_generator
2917       ENDIF
2918       CALL local_stop
2919    ENDIF
2920
2921!
2922!-- Determine damping level index for 1D model
2923    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
2924       IF ( damp_level_1d == -1.0 )  THEN
2925          damp_level_1d     = zu(nzt+1)
2926          damp_level_ind_1d = nzt + 1
2927       ELSEIF ( damp_level_1d < 0.0  .OR.  damp_level_1d > zu(nzt+1) )  THEN
2928          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  damp_level_1d=', &
2929                                    damp_level_1d, ' must be > 0.0 and < ',  &
2930                                    'zu(nzt+1)', zu(nzt+1)
2931          CALL local_stop
2932       ELSE
2933          DO  k = 1, nzt+1
2934             IF ( damp_level_1d <= zu(k) )  THEN
2935                damp_level_ind_1d = k
2936                EXIT
2937             ENDIF
2938          ENDDO
2939       ENDIF
2940    ENDIF
2941!
2942!-- Check some other 1d-model parameters
2943    IF ( TRIM( mixing_length_1d ) /= 'as_in_3d_model'  .AND. &
2944         TRIM( mixing_length_1d ) /= 'blackadar' )  THEN
2945       IF ( myid == 0 )  PRINT*, '+++ check_parameters: mixing_length_1d = "', &
2946                                 TRIM( mixing_length_1d ), '" is unknown'
2947       CALL local_stop
2948    ENDIF
2949    IF ( TRIM( dissipation_1d ) /= 'as_in_3d_model'  .AND. &
2950         TRIM( dissipation_1d ) /= 'detering' )  THEN
2951       IF ( myid == 0 )  PRINT*, '+++ check_parameters: dissipation_1d = "', &
2952                                 TRIM( dissipation_1d ), '" is unknown'
2953       CALL local_stop
2954    ENDIF
2955
2956!
2957!-- Set time for the next user defined restart (time_restart is the
2958!-- internal parameter for steering restart events)
2959    IF ( restart_time /= 9999999.9 )  THEN
2960       IF ( restart_time > simulated_time )  time_restart = restart_time
2961    ELSE
2962!
2963!--    In case of a restart run, set internal parameter to default (no restart)
2964!--    if the NAMELIST-parameter restart_time is omitted
2965       time_restart = 9999999.9
2966    ENDIF
2967
2968!
2969!-- Set default value of the time needed to terminate a model run
2970    IF ( termination_time_needed == -1.0 )  THEN
2971       IF ( host(1:3) == 'ibm' )  THEN
2972          termination_time_needed = 300.0
2973       ELSE
2974          termination_time_needed = 35.0
2975       ENDIF
2976    ENDIF
2977
2978!
2979!-- Check the time needed to terminate a model run
2980    IF ( host(1:3) == 't3e' )  THEN
2981!
2982!--    Time needed must be at least 30 seconds on all CRAY machines, because
2983!--    MPP_TREMAIN gives the remaining CPU time only in steps of 30 seconds
2984       IF ( termination_time_needed <= 30.0 )  THEN
2985          IF ( myid == 0 )  THEN
2986             PRINT*, '+++ check_parameters:  termination_time_needed', &
2987                      termination_time_needed
2988             PRINT*, '                       must be > 30.0 on host "', host, &
2989                     '"'
2990          ENDIF
2991          CALL local_stop
2992       ENDIF
2993    ELSEIF ( host(1:3) == 'ibm' )  THEN
2994!
2995!--    On IBM-regatta machines the time should be at least 300 seconds,
2996!--    because the job time consumed before executing palm (for compiling,
2997!--    copying of files, etc.) has to be regarded
2998       IF ( termination_time_needed < 300.0 )  THEN
2999          IF ( myid == 0 )  THEN
3000             PRINT*, '+++ WARNING: check_parameters:  termination_time_',  &
3001                         'needed', termination_time_needed
3002             PRINT*, '                                should be >= 300.0', &
3003                         ' on host "', host, '"'
3004          ENDIF
3005       ENDIF
3006    ENDIF
3007
3008
3009 END SUBROUTINE check_parameters
Note: See TracBrowser for help on using the repository browser.