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

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

new dvrp features added

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