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

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