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

Last change on this file since 133 was 132, checked in by letzel, 16 years ago

Vertical profiles now based on nzb_s_inner; they are divided by
ngp_2dh_s_inner (scalars, procucts of scalars and velocity components) and
ngp_2dh (staggered velocity components and their products), respectively.

Allow new case bc_uv_t = 'dirichlet_0' for channel flow.

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