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

Last change on this file since 169 was 153, checked in by steinfeld, 17 years ago

Update for the plant canopy model

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