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

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

further preliminary changes concerning coupling

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