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

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

preliminary version, several changes to be explained later

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