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

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