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

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

preliminary changes for radiation conditions

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