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

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

preliminary changes for precipitation output

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