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

Last change on this file since 335 was 332, checked in by raasch, 15 years ago

bugfix in check_parameters concerning ocean version

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