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

Last change on this file since 531 was 531, checked in by heinze, 14 years ago

call of subsidence added in prognostic equations for humidity/passive scalar, some bugfixes

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