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

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