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

Last change on this file since 75 was 75, checked in by raasch, 18 years ago

preliminary update for changes concerning non-cyclic boundary conditions

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