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

Last change on this file since 197 was 197, checked in by raasch, 16 years ago

further adjustments for SGI and other small changes

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