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

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

more preliminary uncomplete changes for ocean version

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