source: palm/tags/release-3.6/SOURCE/check_parameters.f90 @ 4417

Last change on this file since 4417 was 232, checked in by raasch, 15 years ago

bugfix concerning output of p profile

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