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

Last change on this file since 109 was 109, checked in by letzel, 17 years ago
  • Bugfix in surface_coupler
  • mrun: completely remove workaround on lcfimm to propagate environment

variables out to the nodes in coupled mode

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