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

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

updating parts of Marcus changes

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