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

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

preliminary changes concerning update of BC-scheme

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