source: palm/tags/release-3.5/SOURCE/check_parameters.f90 @ 203

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

file headers updated for the next release 3.5

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