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

Last change on this file since 138 was 138, checked in by letzel, 16 years ago

Plant canopy model of Watanabe (2004,BLM 112,307-341) added.

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