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

Last change on this file since 60 was 57, checked in by raasch, 17 years ago

preliminary update of further changes, advec_particles is not running!

  • Property svn:keywords set to Id
File size: 83.2 KB
Line 
1 SUBROUTINE check_parameters
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6! "by_user" allowed as initializing action, -data_output_ts,
7! leapfrog with non-flat topography not allowed any more, pt_reference is
8! checked
9!
10! Former revisions:
11! -----------------
12! $Id: check_parameters.f90 57 2007-03-09 12:05:41Z raasch $
13!
14! 20 2007-02-26 00:12:32Z raasch
15! Temperature and humidity gradients at top are now calculated for nzt+1,
16! top_heatflux and respective boundary condition bc_pt_t is checked
17!
18! RCS Log replace by Id keyword, revision history cleaned up
19!
20! Revision 1.61  2006/08/04 14:20:25  raasch
21! do2d_unit and do3d_unit now defined as 2d-arrays, check of
22! use_upstream_for_tke, default value for dt_dopts,
23! generation of file header moved from routines palm and header to here
24!
25! Revision 1.1  1997/08/26 06:29:23  raasch
26! Initial revision
27!
28!
29! Description:
30! ------------
31! Check control parameters and deduce further quantities.
32!------------------------------------------------------------------------------!
33
34    USE arrays_3d
35    USE constants
36    USE control_parameters
37    USE grid_variables
38    USE indices
39    USE model_1d
40    USE netcdf_control
41    USE particle_attributes
42    USE pegrid
43    USE profil_parameter
44    USE statistics
45    USE transpose_indices
46
47    IMPLICIT NONE
48
49    CHARACTER (LEN=1)   ::  sq
50    CHARACTER (LEN=6)   ::  var
51    CHARACTER (LEN=7)   ::  unit
52    CHARACTER (LEN=8)   ::  date
53    CHARACTER (LEN=10)  ::  time
54    CHARACTER (LEN=100) ::  action
55
56    INTEGER ::  i, ilen, intervals, iter, j, k, nnxh, nnyh, position, prec
57    LOGICAL ::  found, ldum
58    REAL    ::  gradient, maxn, maxp
59
60
61!
62!-- Warning, if host is not set
63    IF ( host(1:1) == ' ' )  THEN
64       IF ( myid == 0 )  THEN
65          PRINT*, '+++ WARNING: check_parameters:'
66          PRINT*, '    "host" is not set. Please check that environment', &
67                       ' variable "localhost"'
68          PRINT*, '    is set before running PALM'
69       ENDIF
70    ENDIF
71
72!
73!-- Generate the file header which is used as a header for most of PALM's
74!-- output files
75    CALL DATE_AND_TIME( date, time )
76    run_date = date(7:8)//'-'//date(5:6)//'-'//date(3:4)
77    run_time = time(1:2)//':'//time(3:4)//':'//time(5:6)
78
79    WRITE ( run_description_header, '(A,2X,A,A,A,I2.2,2X,A,A,2X,A,1X,A)' )  &
80              TRIM( version ),'run: ', TRIM( run_identifier ), '.', &
81              runnr, 'host: ', TRIM( host ), run_date, run_time
82
83!
84!-- Check topography setting (check for illegal parameter combinations)
85    IF ( topography /= 'flat' )  THEN
86       action = ' '
87       IF ( scalar_advec /= 'pw-scheme' )  THEN
88          WRITE( action, '(A,A)' )  'scalar_advec = ', scalar_advec
89       ENDIF
90       IF ( momentum_advec /= 'pw-scheme' )  THEN
91          WRITE( action, '(A,A)' )  'momentum_advec = ', momentum_advec
92       ENDIF
93       IF ( timestep_scheme(1:8) == 'leapfrog' )  THEN
94          WRITE( action, '(A,A)' )  'timestep_scheme = ', timestep_scheme
95       ENDIF
96       IF ( psolver == 'multigrid'  .OR.  psolver == 'sor' )  THEN
97          WRITE( action, '(A,A)' )  'psolver = ', psolver
98       ENDIF
99       IF ( sloping_surface )  THEN
100          WRITE( action, '(A)' )  'sloping surface = .TRUE.'
101       ENDIF
102       IF ( galilei_transformation )  THEN
103          WRITE( action, '(A)' )  'galilei_transformation = .TRUE.'
104       ENDIF
105       IF ( cloud_physics )  THEN
106          WRITE( action, '(A)' )  'cloud_physics = .TRUE.'
107       ENDIF
108       IF ( cloud_droplets )  THEN
109          WRITE( action, '(A)' )  'cloud_droplets = .TRUE.'
110       ENDIF
111       IF ( moisture )  THEN
112          WRITE( action, '(A)' )  'moisture = .TRUE.'
113       ENDIF
114       IF ( .NOT. prandtl_layer )  THEN
115          WRITE( action, '(A)' )  'prandtl_layer = .FALSE.'
116       ENDIF
117       IF ( action /= ' ' )  THEN
118          IF ( myid == 0 )  THEN
119             PRINT*, '+++ check_parameters:'
120             PRINT*, '    a non-flat topography does not allow ', TRIM( action )
121          ENDIF
122          CALL local_stop
123       ENDIF
124    ENDIF
125!
126!-- Check whether there are any illegal values
127!-- Pressure solver:
128    IF ( psolver /= 'poisfft'  .AND.  psolver /= 'poisfft_hybrid'  .AND. &
129         psolver /= 'sor'  .AND.  psolver /= 'multigrid' )  THEN
130       IF ( myid == 0 )  THEN
131          PRINT*, '+++ check_parameters:'
132          PRINT*, '    unknown solver for perturbation pressure: psolver=', &
133                  psolver
134       ENDIF
135       CALL local_stop
136    ENDIF
137
138#if defined( __parallel )
139    IF ( psolver == 'poisfft_hybrid'  .AND.  pdims(2) /= 1 )  THEN
140       IF ( myid == 0 )  THEN
141          PRINT*, '+++ check_parameters:'
142          PRINT*, '    psolver="', TRIM( psolver ), '" only works for a ', &
143                       '1d domain-decomposition along x'
144          PRINT*, '    please do not set npey/=1 in the parameter file'
145       ENDIF
146       CALL local_stop
147    ENDIF
148    IF ( ( psolver == 'poisfft_hybrid'  .OR.  psolver == 'multigrid' )  .AND.  &
149         ( nxra > nxr  .OR.  nyna > nyn  .OR.  nza > nz ) )  THEN
150       IF ( myid == 0 )  THEN
151          PRINT*, '+++ check_parameters:'
152          PRINT*, '    psolver="', TRIM( psolver ), '" does not work for ', &
153                       'subdomains with unequal size'
154          PRINT*, '    please set grid_matching = ''strict'' in the parameter',&
155                       ' file'
156       ENDIF
157       CALL local_stop
158    ENDIF
159#else
160    IF ( psolver == 'poisfft_hybrid' )  THEN
161       IF ( myid == 0 )  THEN
162          PRINT*, '+++ check_parameters:'
163          PRINT*, '    psolver="', TRIM( psolver ), '" only works for a ', &
164                       'parallel environment'
165       ENDIF
166       CALL local_stop
167    ENDIF
168#endif
169
170    IF ( psolver == 'multigrid' )  THEN
171       IF ( cycle_mg == 'w' )  THEN
172          gamma_mg = 2
173       ELSEIF ( cycle_mg == 'v' )  THEN
174          gamma_mg = 1
175       ELSE
176          IF ( myid == 0 )  THEN
177             PRINT*, '+++ check_parameters:'
178             PRINT*, '    unknown multigrid cycle: cycle_mg=', cycle_mg
179          ENDIF
180          CALL local_stop
181       ENDIF
182    ENDIF
183
184    IF ( fft_method /= 'singleton-algorithm'  .AND.  &
185         fft_method /= 'temperton-algorithm'  .AND.  &
186         fft_method /= 'system-specific' )  THEN
187       IF ( myid == 0 )  THEN
188          PRINT*, '+++ check_parameters:'
189          PRINT*, '    unknown fft-algorithm: fft_method=', fft_method
190       ENDIF
191       CALL local_stop
192    ENDIF
193
194!
195!-- Advection schemes:
196    IF ( momentum_advec /= 'pw-scheme' .AND. momentum_advec /= 'ups-scheme' ) &
197    THEN
198       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  unknown advection ', &
199                                 'scheme: momentum_advec=', momentum_advec
200       CALL local_stop
201    ENDIF
202    IF ( ( momentum_advec == 'ups-scheme'  .OR.  scalar_advec == 'ups-scheme' )&
203                                      .AND.  timestep_scheme /= 'euler' )  THEN
204       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  momentum_advec=', &
205                                 momentum_advec, ' is not allowed with ', &
206                                 'timestep_scheme=', timestep_scheme
207       CALL local_stop
208    ENDIF
209
210    IF ( scalar_advec /= 'pw-scheme'  .AND.  scalar_advec /= 'bc-scheme'  .AND.&
211         scalar_advec /= 'ups-scheme' )  THEN
212       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  unknown advection ', &
213                                 'scheme: scalar_advec=', scalar_advec
214       CALL local_stop
215    ENDIF
216
217    IF ( use_sgs_for_particles  .AND.  .NOT. use_upstream_for_tke )  THEN
218       use_upstream_for_tke = .TRUE.
219       IF ( myid == 0 )  THEN
220          PRINT*, '+++ WARNING: check_parameters:  use_upstream_for_tke set ', &
221                       '.TRUE. because use_sgs_for_particles = .TRUE.'
222       ENDIF
223    ENDIF
224
225    IF ( use_upstream_for_tke  .AND.  timestep_scheme(1:8) == 'leapfrog' )  THEN
226       IF ( myid == 0 )  THEN
227          PRINT*, '+++ check_parameters:  use_upstream_for_tke = .TRUE. ', &
228                       'not allowed with timestep_scheme=', timestep_scheme
229       ENDIF
230       CALL local_stop
231    ENDIF
232
233!
234!-- Timestep schemes:
235    SELECT CASE ( TRIM( timestep_scheme ) )
236
237       CASE ( 'euler' )
238          intermediate_timestep_count_max = 1
239          asselin_filter_factor           = 0.0
240
241       CASE ( 'leapfrog', 'leapfrog+euler' )
242          intermediate_timestep_count_max = 1
243
244       CASE ( 'runge-kutta-2' )
245          intermediate_timestep_count_max = 2
246          asselin_filter_factor           = 0.0
247
248       CASE ( 'runge-kutta-3' )
249          intermediate_timestep_count_max = 3
250          asselin_filter_factor           = 0.0
251
252       CASE DEFAULT
253          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  unknown timestep ',&
254                                    'scheme: timestep_scheme=', timestep_scheme
255          CALL local_stop
256
257    END SELECT
258
259    IF ( scalar_advec /= 'pw-scheme'  .AND.  timestep_scheme(1:5) == 'runge' ) &
260    THEN
261       IF ( myid == 0 )  THEN
262          PRINT*, '+++ check_parameters:  scalar advection scheme "', &
263                                          TRIM( scalar_advec ), '"'
264          PRINT*, '    does not work with timestep_scheme "', &
265                                          TRIM( timestep_scheme ), '"'
266       ENDIF
267       CALL local_stop
268    ENDIF
269
270    IF ( momentum_advec /= 'pw-scheme' .AND. timestep_scheme(1:5) == 'runge' ) &
271    THEN
272       IF ( myid == 0 )  THEN
273          PRINT*, '+++ check_parameters:  momentum advection scheme "', &
274                                          TRIM( momentum_advec ), '"'
275          PRINT*, '    does not work with timestep_scheme "', &
276                                          TRIM( timestep_scheme ), '"'
277       ENDIF
278       CALL local_stop
279    ENDIF
280
281
282    IF ( initializing_actions == ' ' )  THEN
283       IF ( myid == 0 )  THEN
284          PRINT*, '+++ check parameters:'
285          PRINT*, '    no value found for initializing_actions'
286       ENDIF
287       CALL local_stop
288    ENDIF
289
290    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
291!
292!--    No model continuation run; several initialising actions are possible
293       action = initializing_actions
294       DO WHILE ( TRIM( action ) /= '' )
295          position = INDEX( action, ' ' )
296          SELECT CASE ( action(1:position-1) )
297
298             CASE ( 'set_constant_profiles', 'set_1d-model_profiles', &
299                    'by_user', 'initialize_vortex',     'initialize_ptanom' )
300                action = action(position+1:)
301
302             CASE DEFAULT
303                IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializi', &
304                               'ng_action unkown or not allowed: action = "', &
305                               TRIM(action), '"'
306                CALL local_stop
307
308          END SELECT
309       ENDDO
310    ENDIF
311    IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0  .AND. &
312         INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
313       IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializing_actions', &
314          '"set_constant_profiles" and "set_1d-model_profiles" are not', &
315          ' allowed simultaneously'
316       CALL local_stop
317    ENDIF
318    IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0  .AND. &
319         INDEX( initializing_actions, 'by_user' ) /= 0 )  THEN
320       IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializing_actions', &
321          '"set_constant_profiles" and "by_user" are not', &
322          ' allowed simultaneously'
323       CALL local_stop
324    ENDIF
325    IF ( INDEX( initializing_actions, 'by_user' ) /= 0  .AND. &
326         INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
327       IF ( myid == 0 )  PRINT*, '+++ check_parameters: initializing_actions', &
328          '"by_user" and "set_1d-model_profiles" are not', &
329          ' allowed simultaneously'
330       CALL local_stop
331    ENDIF
332
333    IF ( cloud_physics  .AND.  .NOT. moisture )  THEN
334       IF ( myid == 0 )  PRINT*, '+++ check_parameters: moisture =', &
335                                 moisture, ' is not allowed with ',  &
336                                 'cloud_physics=', cloud_physics
337       CALL local_stop
338    ENDIF
339
340    IF ( moisture  .AND.  sloping_surface )  THEN
341       IF ( myid == 0 )  PRINT*, '+++ check_parameters: moisture = TRUE', &
342                                 'and hang = TRUE are not',               &
343                                 ' allowed simultaneously' 
344       CALL local_stop       
345    ENDIF
346
347    IF ( moisture  .AND.  scalar_advec == 'ups-scheme' )  THEN
348       IF ( myid == 0 )  PRINT*, '+++ check_parameters: UPS-scheme', &
349                                 'is not implemented for moisture'
350       CALL local_stop       
351    ENDIF
352
353    IF ( passive_scalar  .AND.  moisture )  THEN
354       IF ( myid == 0 )  PRINT*, '+++ check_parameters: moisture = TRUE and', &
355                                 'passive_scalar = TRUE is not allowed ',     &
356                                 'simultaneously'
357       CALL local_stop
358    ENDIF
359
360    IF ( passive_scalar  .AND.  scalar_advec == 'ups-scheme' )  THEN
361       IF ( myid == 0 )  PRINT*, '+++ check_parameters: UPS-scheme', &
362                                 'is not implemented for passive_scalar'
363       CALL local_stop       
364    ENDIF
365
366    IF ( grid_matching /= 'strict'  .AND.  grid_matching /= 'match' )  THEN
367       IF ( myid == 0 )  PRINT*, '+++ check_parameters: illegal value "', &
368                                 TRIM( grid_matching ),                   &
369                                 '" found for parameter grid_matching'
370       CALL local_stop       
371    ENDIF
372
373!
374!-- In case of no model continuation run, check initialising parameters and
375!-- deduce further quantities
376    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
377
378!
379!--    Initial profiles for 1D and 3D model, respectively
380       u_init  = ug_surface
381       v_init  = vg_surface
382       pt_init = pt_surface
383       IF ( moisture )        q_init = q_surface
384       IF ( passive_scalar )  q_init = s_surface
385
386!
387!--
388!--    If required, compute initial profile of the geostrophic wind
389!--    (component ug)
390       i = 1
391       gradient = 0.0
392       ug_vertical_gradient_level_ind(1) = 0
393       ug(0) = ug_surface
394       DO  k = 1, nzt+1
395          IF ( ug_vertical_gradient_level(i) < zu(k)  .AND. &
396               ug_vertical_gradient_level(i) >= 0.0 )  THEN
397             gradient = ug_vertical_gradient(i) / 100.0
398             ug_vertical_gradient_level_ind(i) = k - 1
399             i = i + 1
400             IF ( i > 10 )  THEN
401                IF ( myid == 0 )  THEN
402                   PRINT*, '+++ check_parameters: upper bound 10 of array', &
403                           ' "ug_vertical_gradient_level_ind" exceeded'
404                ENDIF
405                CALL local_stop
406             ENDIF
407          ENDIF
408          IF ( gradient /= 0.0 )  THEN
409             IF ( k /= 1 )  THEN
410                ug(k) = ug(k-1) + dzu(k) * gradient
411             ELSE
412                ug(k) = ug_surface + 0.5 * dzu(k) * gradient
413             ENDIF
414          ELSE
415             ug(k) = ug(k-1)
416          ENDIF
417       ENDDO
418
419       u_init = ug
420
421!
422!--    In case of no given gradients for ug, choose a vanishing gradient
423       IF ( ug_vertical_gradient_level(1) == -1.0 )  THEN
424          ug_vertical_gradient_level(1) = 0.0
425       ENDIF 
426
427!
428!--
429!--    If required, compute initial profile of the geostrophic wind
430!--    (component vg)
431       i = 1
432       gradient = 0.0
433       vg_vertical_gradient_level_ind(1) = 0
434       vg(0) = vg_surface
435       DO  k = 1, nzt+1
436          IF ( vg_vertical_gradient_level(i) < zu(k)  .AND. &
437               vg_vertical_gradient_level(i) >= 0.0 )  THEN
438             gradient = vg_vertical_gradient(i) / 100.0
439             vg_vertical_gradient_level_ind(i) = k - 1
440             i = i + 1
441             IF ( i > 10 )  THEN
442                IF ( myid == 0 )  THEN
443                   PRINT*, '+++ check_parameters: upper bound 10 of array', &
444                           ' "vg_vertical_gradient_level_ind" exceeded'
445                ENDIF
446                CALL local_stop
447             ENDIF
448          ENDIF
449          IF ( gradient /= 0.0 )  THEN
450             IF ( k /= 1 )  THEN
451                vg(k) = vg(k-1) + dzu(k) * gradient
452             ELSE
453                vg(k) = vg_surface + 0.5 * dzu(k) * gradient
454             ENDIF
455          ELSE
456             vg(k) = vg(k-1)
457          ENDIF
458       ENDDO
459
460       v_init = vg
461 
462!
463!--    In case of no given gradients for vg, choose a vanishing gradient
464       IF ( vg_vertical_gradient_level(1) == -1.0 )  THEN
465          vg_vertical_gradient_level(1) = 0.0
466       ENDIF
467
468!
469!--    If required, compute initial temperature profile using the given
470!--    temperature gradient
471       i = 1
472       gradient = 0.0
473       pt_vertical_gradient_level_ind(1) = 0
474       DO  k = 1, nzt+1
475          IF ( pt_vertical_gradient_level(i) < zu(k)  .AND. &
476               pt_vertical_gradient_level(i) >= 0.0 )  THEN
477             gradient = pt_vertical_gradient(i) / 100.0
478             pt_vertical_gradient_level_ind(i) = k - 1
479             i = i + 1
480             IF ( i > 10 )  THEN
481                IF ( myid == 0 )  THEN
482                   PRINT*, '+++ check_parameters: upper bound 10 of array', &
483                           ' "pt_vertical_gradient_level_ind" exceeded'
484                ENDIF
485                CALL local_stop
486             ENDIF
487          ENDIF
488          IF ( gradient /= 0.0 )  THEN
489             IF ( k /= 1 )  THEN
490                pt_init(k) = pt_init(k-1) + dzu(k) * gradient
491             ELSE
492                pt_init(k) = pt_init(k-1) + 0.5 * dzu(k) * gradient
493             ENDIF
494          ELSE
495             pt_init(k) = pt_init(k-1)
496          ENDIF
497       ENDDO
498
499!
500!--    In case of no given temperature gradients, choose gradient of neutral
501!--    stratification
502       IF ( pt_vertical_gradient_level(1) == -1.0 )  THEN
503          pt_vertical_gradient_level(1) = 0.0
504       ENDIF
505
506!
507!--    Store temperature gradient at the top boundary for possile Neumann
508!--    boundary condition
509       bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1)
510
511!
512!--    If required, compute initial humidity or scalar profile using the given
513!--    humidity/scalar gradient. In case of scalar transport, initially store
514!--    values of the scalar parameters on humidity parameters
515       IF ( passive_scalar )  THEN
516          bc_q_b                    = bc_s_b
517          bc_q_t                    = bc_s_t
518          q_surface                 = s_surface
519          q_surface_initial_change  = s_surface_initial_change
520          q_vertical_gradient       = s_vertical_gradient
521          q_vertical_gradient_level = s_vertical_gradient_level
522          surface_waterflux         = surface_scalarflux
523       ENDIF
524
525       IF ( moisture  .OR.  passive_scalar )  THEN
526
527          i = 1
528          gradient = 0.0
529          q_vertical_gradient_level_ind(1) = 0
530          DO  k = 1, nzt+1
531             IF ( q_vertical_gradient_level(i) < zu(k)  .AND. &
532                  q_vertical_gradient_level(i) >= 0.0 )  THEN
533                gradient = q_vertical_gradient(i) / 100.0
534                q_vertical_gradient_level_ind(i) = k - 1
535                i = i + 1
536                IF ( i > 10 )  THEN
537                   IF ( myid == 0 )  THEN
538                      PRINT*, '+++ check_parameters: upper bound 10 of arr', &
539                              'ay "q_vertical_gradient_level_ind" exceeded'
540                   ENDIF
541                   CALL local_stop
542                ENDIF
543             ENDIF
544             IF ( gradient /= 0.0 )  THEN
545                IF ( k /= 1 )  THEN
546                   q_init(k) = q_init(k-1) + dzu(k) * gradient
547                ELSE
548                   q_init(k) = q_init(k-1) + 0.5 * dzu(k) * gradient
549                ENDIF
550             ELSE
551                q_init(k) = q_init(k-1)
552             ENDIF
553          ENDDO
554
555!
556!--       In case of no given humidity gradients, choose zero gradient
557!--       conditions
558          IF ( q_vertical_gradient_level(1) == -1.0 )  THEN
559             q_vertical_gradient_level(1) = 0.0
560          ENDIF
561
562!
563!--       Store humidity gradient at the top boundary for possile Neumann
564!--       boundary condition
565          bc_q_t_val = ( q_init(nzt+1) - q_init(nzt) ) / dzu(nzt+1)
566
567       ENDIF
568
569    ENDIF
570
571!
572!-- Compute Coriolis parameter
573    f  = 2.0 * omega * SIN( phi / 180.0 * pi )
574    fs = 2.0 * omega * COS( phi / 180.0 * pi )
575
576!
577!-- Reference temperature to be used in buoyancy terms
578    IF ( pt_reference /= 9999999.9 )  use_pt_reference = .TRUE.
579
580!
581!-- In case of a given slope, compute the relevant quantities
582    IF ( alpha_surface /= 0.0 )  THEN
583       IF ( ABS( alpha_surface ) > 90.0 )  THEN
584          IF ( myid == 0 )  PRINT*, '+++ check_parameters: ABS( alpha_surface',&
585                                    '=', alpha_surface, ' ) must be < 90.0'
586          CALL local_stop
587       ENDIF
588       sloping_surface = .TRUE.
589       cos_alpha_surface = COS( alpha_surface / 180.0 * pi )
590       sin_alpha_surface = SIN( alpha_surface / 180.0 * pi )
591    ENDIF
592
593!
594!-- Check time step and cfl_factor
595    IF ( dt /= -1.0 )  THEN
596       IF ( dt <= 0.0  .AND.  dt /= -1.0 )  THEN
597          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  dt=', dt, ' <= 0.0'
598          CALL local_stop
599       ENDIF
600       dt_3d = dt
601       dt_fixed = .TRUE.
602    ENDIF
603
604    IF ( cfl_factor <= 0.0  .OR.  cfl_factor > 1.0 )  THEN
605       IF ( cfl_factor == -1.0 )  THEN
606          IF ( momentum_advec == 'ups-scheme'  .OR.  &
607               scalar_advec == 'ups-scheme' )  THEN
608             cfl_factor = 0.8
609          ELSE
610             IF ( timestep_scheme == 'runge-kutta-2' )  THEN
611                cfl_factor = 0.8
612             ELSEIF ( timestep_scheme == 'runge-kutta-3' )  THEN
613                cfl_factor = 0.9
614             ELSE
615                cfl_factor = 0.1
616             ENDIF
617          ENDIF
618       ELSE
619          IF ( myid == 0 )  THEN
620             PRINT*, '+++ check_parameters: cfl_factor=', cfl_factor, &
621                         ' out of range'
622             PRINT*, '+++                   0.0 < cfl_factor <= 1.0 is required'
623          ENDIF
624          CALL local_stop
625       ENDIF
626    ENDIF
627
628!
629!-- Store simulated time at begin
630    simulated_time_at_begin = simulated_time
631
632!
633!-- Set wind speed in the Galilei-transformed system
634    IF ( galilei_transformation )  THEN
635       IF ( use_ug_for_galilei_tr .AND.                &
636            ug_vertical_gradient_level(1) == 0.0 .AND. & 
637            vg_vertical_gradient_level(1) == 0.0 )  THEN
638          u_gtrans = ug_surface
639          v_gtrans = vg_surface
640       ELSEIF ( use_ug_for_galilei_tr .AND.                &
641                ug_vertical_gradient_level(1) /= 0.0 )  THEN
642          IF ( myid == 0 )  THEN
643             PRINT*, '+++ check_parameters:'
644             PRINT*, '    baroclinicity (ug) not allowed'
645             PRINT*, '    simultaneously with galilei transformation'
646          ENDIF
647          CALL local_stop
648       ELSEIF ( use_ug_for_galilei_tr .AND.                &
649                vg_vertical_gradient_level(1) /= 0.0 )  THEN
650          IF ( myid == 0 )  THEN
651             PRINT*, '+++ check_parameters:'
652             PRINT*, '    baroclinicity (vg) not allowed'
653             PRINT*, '    simultaneously with galilei transformation'
654          ENDIF
655          CALL local_stop
656       ELSE
657          IF ( myid == 0 )  THEN
658             PRINT*, '+++ WARNING: check_parameters:'
659             PRINT*, '    variable translation speed used for galilei-tran' // &
660                          'sformation, which'
661             PRINT*, '    may cause instabilities in stably stratified regions'
662          ENDIF
663       ENDIF
664    ENDIF
665
666!
667!-- In case of using a prandtl-layer, calculated (or prescribed) surface
668!-- fluxes have to be used in the diffusion-terms
669    IF ( prandtl_layer )  use_surface_fluxes = .TRUE.
670
671!
672!-- Check boundary conditions and set internal variables:
673!-- Lateral boundary conditions
674    IF ( bc_lr /= 'cyclic'  .AND.  bc_lr /= 'dirichlet/neumann'  .AND. &
675         bc_lr /= 'neumann/dirichlet' )  THEN
676       IF ( myid == 0 )  THEN
677          PRINT*, '+++ check_parameters:'
678          PRINT*, '    unknown boundary condition: bc_lr = ', bc_lr
679       ENDIF
680       CALL local_stop
681    ENDIF
682    IF ( bc_ns /= 'cyclic'  .AND.  bc_ns /= 'dirichlet/neumann'  .AND. &
683         bc_ns /= 'neumann/dirichlet' )  THEN
684       IF ( myid == 0 )  THEN
685          PRINT*, '+++ check_parameters:'
686          PRINT*, '    unknown boundary condition: bc_ns = ', bc_ns
687       ENDIF
688       CALL local_stop
689    ENDIF
690
691!
692!-- Non-cyclic lateral boundaries require the multigrid method and Piascek-
693!-- Willimas advection scheme. Several schemes and tools do not work with
694!-- non-cyclic boundary conditions.
695    IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
696       IF ( psolver /= 'multigrid' )  THEN
697          IF ( myid == 0 )  THEN
698             PRINT*, '+++ check_parameters:'
699             PRINT*, '    non-cyclic lateral boundaries do not allow', &
700                          ' psolver = ', psolver 
701          ENDIF
702          CALL local_stop
703       ENDIF
704       IF ( momentum_advec /= 'pw-scheme' )  THEN
705          IF ( myid == 0 )  THEN
706             PRINT*, '+++ check_parameters:'
707             PRINT*, '    non-cyclic lateral boundaries do not allow', &
708                          ' momentum_advec = ', momentum_advec 
709          ENDIF
710          CALL local_stop
711       ENDIF
712       IF ( scalar_advec /= 'pw-scheme' )  THEN
713          IF ( myid == 0 )  THEN
714             PRINT*, '+++ check_parameters:'
715             PRINT*, '    non-cyclic lateral boundaries do not allow', &
716                          ' scalar_advec = ', scalar_advec 
717          ENDIF
718          CALL local_stop
719       ENDIF
720       IF ( galilei_transformation )  THEN
721          IF ( myid == 0 )  THEN
722             PRINT*, '+++ check_parameters:'
723             PRINT*, '    non-cyclic lateral boundaries do not allow', &
724                          ' galilei_transformation = .T.' 
725          ENDIF
726          CALL local_stop
727       ENDIF
728       IF ( conserve_volume_flow )  THEN
729          IF ( myid == 0 )  THEN
730             PRINT*, '+++ check_parameters:'
731             PRINT*, '    non-cyclic lateral boundaries do not allow', &
732                          ' conserve_volume_flow = .T.' 
733          ENDIF
734          CALL local_stop
735       ENDIF
736    ENDIF
737
738!
739!-- Bottom boundary condition for the turbulent Kinetic energy
740    IF ( bc_e_b == 'neumann' )  THEN
741       ibc_e_b = 1
742       IF ( adjust_mixing_length  .AND.  prandtl_layer )  THEN
743          IF ( myid == 0 )  THEN
744             PRINT*, '+++ WARNING: check_parameters:'
745             PRINT*, '    adjust_mixing_length = TRUE and bc_e_b = ', bc_e_b
746          ENDIF
747       ENDIF
748    ELSEIF ( bc_e_b == '(u*)**2+neumann' )  THEN
749       ibc_e_b = 2
750       IF ( .NOT. adjust_mixing_length  .AND.  prandtl_layer )  THEN
751          IF ( myid == 0 )  THEN
752             PRINT*, '+++ WARNING: check_parameters:'
753             PRINT*, '    adjust_mixing_length = FALSE and bc_e_b = ', bc_e_b
754          ENDIF
755       ENDIF
756       IF ( .NOT. prandtl_layer )  THEN
757          bc_e_b = 'neumann'
758          ibc_e_b = 1
759          IF ( myid == 0 )  THEN
760             PRINT*, '+++ WARNING: check_parameters:'
761             PRINT*, '    boundary condition bc_e_b changed to "', bc_e_b, '"'
762          ENDIF
763       ENDIF
764    ELSE
765       IF ( myid == 0 )  THEN
766          PRINT*, '+++ check_parameters:'
767          PRINT*, '    unknown boundary condition: bc_e_b = ', bc_e_b
768       ENDIF
769       CALL local_stop
770    ENDIF
771
772!
773!-- Boundary conditions for perturbation pressure
774    IF ( bc_p_b == 'dirichlet' )  THEN
775       ibc_p_b = 0
776    ELSEIF ( bc_p_b == 'neumann' )  THEN
777       ibc_p_b = 1
778    ELSEIF ( bc_p_b == 'neumann+inhomo' )  THEN
779       ibc_p_b = 2
780    ELSE
781       IF ( myid == 0 )  THEN
782          PRINT*, '+++ check_parameters:'
783          PRINT*, '    unknown boundary condition: bc_p_b = ', bc_p_b
784       ENDIF
785       CALL local_stop
786    ENDIF
787    IF ( ibc_p_b == 2  .AND.  .NOT. prandtl_layer )  THEN
788       IF ( myid == 0 )  THEN
789          PRINT*, '+++ check_parameters:'
790          PRINT*, '    boundary condition: bc_p_b = ', TRIM( bc_p_b ), &
791                       ' not allowed with'
792          PRINT*, '    prandtl_layer = .FALSE.' 
793       ENDIF
794       CALL local_stop
795    ENDIF
796    IF ( bc_p_t == 'dirichlet' )  THEN
797       ibc_p_t = 0
798    ELSEIF ( bc_p_t == 'neumann' )  THEN
799       ibc_p_t = 1
800    ELSE
801       IF ( myid == 0 )  THEN
802          PRINT*, '+++ check_parameters:'
803          PRINT*, '    unknown boundary condition: bc_p_t = ', bc_p_t
804       ENDIF
805       CALL local_stop
806    ENDIF
807
808!
809!-- Boundary conditions for potential temperature
810    IF ( bc_pt_b == 'dirichlet' )  THEN
811       ibc_pt_b = 0
812    ELSEIF ( bc_pt_b == 'neumann' )  THEN
813       ibc_pt_b = 1
814    ELSE
815       IF ( myid == 0 )  THEN
816          PRINT*, '+++ check_parameters:'
817          PRINT*, '    unknown boundary condition: bc_pt_b = ', bc_pt_b
818       ENDIF
819       CALL local_stop
820    ENDIF
821    IF ( bc_pt_t == 'dirichlet' )  THEN
822       ibc_pt_t = 0
823    ELSEIF ( bc_pt_t == 'neumann' )  THEN
824       ibc_pt_t = 1
825    ELSEIF ( bc_pt_t == 'initial_gradient' )  THEN
826       ibc_pt_t = 2
827    ELSE
828       IF ( myid == 0 )  THEN
829          PRINT*, '+++ check_parameters:'
830          PRINT*, '    unknown boundary condition: bc_pt_t = ', bc_pt_t
831       ENDIF
832       CALL local_stop
833    ENDIF
834
835    IF ( surface_heatflux == 9999999.9 )  constant_heatflux     = .FALSE.
836    IF ( top_heatflux     == 9999999.9 )  constant_top_heatflux = .FALSE.
837
838!
839!-- A given surface temperature implies Dirichlet boundary condition for
840!-- temperature. In this case specification of a constant heat flux is
841!-- forbidden.
842    IF ( ibc_pt_b == 0  .AND.   constant_heatflux  .AND. &
843         surface_heatflux /= 0.0 )  THEN
844       IF ( myid == 0 )  THEN
845          PRINT*, '+++ check_parameters:'
846          PRINT*, '    boundary_condition: bc_pt_b = ', bc_pt_b
847          PRINT*, '    is not allowed with constant_heatflux = .TRUE.'
848       ENDIF
849       CALL local_stop
850    ENDIF
851    IF ( constant_heatflux  .AND.  pt_surface_initial_change /= 0.0 )  THEN
852       IF ( myid == 0 )  THEN
853          PRINT*, '+++ check_parameters: constant_heatflux = .TRUE. is not'
854          PRINT*, '    allowed with pt_surface_initial_change (/=0) = ', &
855                  pt_surface_initial_change
856       ENDIF
857       CALL local_stop
858    ENDIF
859
860!
861!-- A given temperature at the top implies Dirichlet boundary condition for
862!-- temperature. In this case specification of a constant heat flux is
863!-- forbidden.
864    IF ( ibc_pt_t == 0  .AND.   constant_top_heatflux  .AND. &
865         top_heatflux /= 0.0 )  THEN
866       IF ( myid == 0 )  THEN
867          PRINT*, '+++ check_parameters:'
868          PRINT*, '    boundary_condition: bc_pt_t = ', bc_pt_t
869          PRINT*, '    is not allowed with constant_top_heatflux = .TRUE.'
870       ENDIF
871       CALL local_stop
872    ENDIF
873
874!
875!-- In case of moisture or passive scalar, set boundary conditions for total
876!-- water content / scalar
877    IF ( moisture  .OR.  passive_scalar ) THEN
878       IF ( moisture )  THEN
879          sq = 'q'
880       ELSE
881          sq = 's'
882       ENDIF
883       IF ( bc_q_b == 'dirichlet' )  THEN
884          ibc_q_b = 0
885       ELSEIF ( bc_q_b == 'neumann' )  THEN
886          ibc_q_b = 1
887       ELSE
888          IF ( myid == 0 )  THEN
889             PRINT*, '+++ check_parameters:'
890             PRINT*, '    unknown boundary condition: bc_', sq, '_b = ', bc_q_b
891          ENDIF
892          CALL local_stop
893       ENDIF
894       IF ( bc_q_t == 'dirichlet' )  THEN
895          ibc_q_t = 0
896       ELSEIF ( bc_q_t == 'neumann' )  THEN
897          ibc_q_t = 1
898       ELSE
899          IF ( myid == 0 )  THEN
900             PRINT*, '+++ check_parameters:'
901             PRINT*, '    unknown boundary condition: bc_', sq, '_t = ', bc_q_t
902          ENDIF
903          CALL local_stop
904       ENDIF
905
906       IF ( surface_waterflux == 0.0 )  constant_waterflux = .FALSE.
907
908!
909!--    A given surface humidity implies Dirichlet boundary condition for
910!--    moisture. In this case specification of a constant water flux is
911!--    forbidden.
912       IF ( ibc_q_b == 0  .AND.  constant_waterflux )  THEN
913          IF ( myid == 0 )  THEN
914             PRINT*, '+++ check_parameters:'
915             PRINT*, '    boundary_condition: bc_', sq, '_b = ', bc_q_b
916             PRINT*, '    is not allowed with prescribed surface flux'
917          ENDIF
918          CALL local_stop
919       ENDIF
920       IF ( constant_waterflux  .AND.  q_surface_initial_change /= 0.0 )  THEN
921          IF ( myid == 0 )  THEN
922             PRINT*, '+++ check_parameters: a prescribed surface flux is not'
923             PRINT*, '    allowed with ', sq, '_surface_initial_change (/=0)', &
924                     ' = ', q_surface_initial_change
925          ENDIF
926          CALL local_stop
927       ENDIF
928       
929    ENDIF
930
931!
932!-- Boundary conditions for horizontal components of wind speed
933    IF ( bc_uv_b == 'dirichlet' )  THEN
934       ibc_uv_b = 0
935    ELSEIF ( bc_uv_b == 'neumann' )  THEN
936       ibc_uv_b = 1
937       IF ( prandtl_layer )  THEN
938          IF ( myid == 0 )  THEN
939             PRINT*, '+++ check_parameters:'
940             PRINT*, '    boundary condition: bc_uv_b = ', TRIM( bc_uv_b ), &
941                          ' is not allowed with'
942             PRINT*, '    prandtl_layer = .TRUE.' 
943          ENDIF
944          CALL local_stop
945       ENDIF
946    ELSE
947       IF ( myid == 0 )  THEN
948          PRINT*, '+++ check_parameters:'
949          PRINT*, '    unknown boundary condition: bc_uv_b = ', bc_uv_b
950       ENDIF
951       CALL local_stop
952    ENDIF
953    IF ( bc_uv_t == 'dirichlet' )  THEN
954       ibc_uv_t = 0
955    ELSEIF ( bc_uv_t == 'neumann' )  THEN
956       ibc_uv_t = 1
957    ELSE
958       IF ( myid == 0 )  THEN
959          PRINT*, '+++ check_parameters:'
960          PRINT*, '    unknown boundary condition: bc_uv_t = ', bc_uv_t
961       ENDIF
962       CALL local_stop
963    ENDIF
964
965!
966!-- Compute and check, respectively, the Rayleigh Damping parameter
967    IF ( rayleigh_damping_factor == -1.0 )  THEN
968       IF ( momentum_advec == 'ups-scheme' )  THEN
969          rayleigh_damping_factor = 0.01
970       ELSE
971          rayleigh_damping_factor = 0.0
972       ENDIF
973    ELSE
974       IF ( rayleigh_damping_factor < 0.0 .OR. rayleigh_damping_factor > 1.0 ) &
975       THEN
976          IF ( myid == 0 )  THEN
977             PRINT*, '+++ check_parameters:'
978             PRINT*, '    rayleigh_damping_factor = ', rayleigh_damping_factor,&
979                          ' out of range [0.0,1.0]'
980          ENDIF
981          CALL local_stop
982       ENDIF
983    ENDIF
984
985    IF ( rayleigh_damping_height == -1.0 )  THEN
986       rayleigh_damping_height = 0.66666666666 * zu(nzt)
987    ELSE
988       IF ( rayleigh_damping_height < 0.0  .OR. &
989            rayleigh_damping_height > zu(nzt) )  THEN
990          IF ( myid == 0 )  THEN
991             PRINT*, '+++ check_parameters:'
992             PRINT*, '    rayleigh_damping_height = ', rayleigh_damping_height,&
993                          ' out of range [0.0,', zu(nzt), ']'
994          ENDIF
995          CALL local_stop
996       ENDIF
997    ENDIF
998
999!
1000!-- Check limiters for Upstream-Spline scheme
1001    IF ( overshoot_limit_u < 0.0  .OR.  overshoot_limit_v < 0.0  .OR.  &
1002         overshoot_limit_w < 0.0  .OR.  overshoot_limit_pt < 0.0  .OR. &
1003         overshoot_limit_e < 0.0 )  THEN
1004       IF ( myid == 0 )  THEN
1005          PRINT*, '+++ check_parameters:'
1006          PRINT*, '    overshoot_limit_... < 0.0 is not allowed'
1007       ENDIF
1008       CALL local_stop
1009    ENDIF
1010    IF ( ups_limit_u < 0.0 .OR. ups_limit_v < 0.0 .OR. ups_limit_w < 0.0 .OR. &
1011         ups_limit_pt < 0.0 .OR. ups_limit_e < 0.0 )  THEN
1012       IF ( myid == 0 )  THEN
1013          PRINT*, '+++ check_parameters:'
1014          PRINT*, '    ups_limit_... < 0.0 is not allowed'
1015       ENDIF
1016       CALL local_stop
1017    ENDIF
1018
1019!
1020!-- Check number of chosen statistic regions. More than 10 regions are not
1021!-- allowed, because so far no more than 10 corresponding output files can
1022!-- be opened (cf. check_open)
1023    IF ( statistic_regions > 9  .OR.  statistic_regions < 0 )  THEN
1024       IF ( myid == 0 )  THEN
1025          PRINT*, '+++ check_parameters: Number of statistic_regions = ', &
1026                       statistic_regions+1
1027          PRINT*, '    Only 10 regions are allowed'
1028       ENDIF
1029       CALL local_stop
1030    ENDIF
1031    IF ( normalizing_region > statistic_regions  .OR. &
1032         normalizing_region < 0)  THEN
1033       IF ( myid == 0 )  THEN
1034          PRINT*, '+++ check_parameters: normalizing_region = ', &
1035                       normalizing_region, ' is unknown'
1036          PRINT*, '    Must be <= ', statistic_regions
1037       ENDIF
1038       CALL local_stop
1039    ENDIF
1040
1041!
1042!-- Set the default intervals for data output, if necessary
1043!-- NOTE: dt_dosp has already been set in package_parin
1044    IF ( dt_data_output /= 9999999.9 )  THEN
1045       IF ( dt_dopr           == 9999999.9 )  dt_dopr           = dt_data_output
1046       IF ( dt_dopts          == 9999999.9 )  dt_dopts          = dt_data_output
1047       IF ( dt_do2d_xy        == 9999999.9 )  dt_do2d_xy        = dt_data_output
1048       IF ( dt_do2d_xz        == 9999999.9 )  dt_do2d_xz        = dt_data_output
1049       IF ( dt_do2d_yz        == 9999999.9 )  dt_do2d_yz        = dt_data_output
1050       IF ( dt_do3d           == 9999999.9 )  dt_do3d           = dt_data_output
1051       IF ( dt_data_output_av == 9999999.9 )  dt_data_output_av = dt_data_output
1052    ENDIF
1053
1054!
1055!-- Set the default skip time intervals for data output, if necessary
1056    IF ( skip_time_dopr    == 9999999.9 ) &
1057                                       skip_time_dopr    = skip_time_data_output
1058    IF ( skip_time_dosp    == 9999999.9 ) &
1059                                       skip_time_dosp    = skip_time_data_output
1060    IF ( skip_time_do2d_xy == 9999999.9 ) &
1061                                       skip_time_do2d_xy = skip_time_data_output
1062    IF ( skip_time_do2d_xz == 9999999.9 ) &
1063                                       skip_time_do2d_xz = skip_time_data_output
1064    IF ( skip_time_do2d_yz == 9999999.9 ) &
1065                                       skip_time_do2d_yz = skip_time_data_output
1066    IF ( skip_time_do3d    == 9999999.9 ) &
1067                                       skip_time_do3d    = skip_time_data_output
1068    IF ( skip_time_data_output_av == 9999999.9 ) &
1069                                skip_time_data_output_av = skip_time_data_output
1070
1071!
1072!-- Check the average intervals (first for 3d-data, then for profiles and
1073!-- spectra)
1074    IF ( averaging_interval > dt_data_output_av )  THEN
1075       IF ( myid == 0 )  THEN
1076          PRINT*, '+++ check_parameters: average_interval=',              &
1077                       averaging_interval, ' must be <= dt_data_output=', &
1078                       dt_data_output
1079       ENDIF
1080       CALL local_stop
1081    ENDIF
1082
1083    IF ( averaging_interval_pr == 9999999.9 )  THEN
1084       averaging_interval_pr = averaging_interval
1085    ENDIF
1086
1087    IF ( averaging_interval_pr > dt_dopr )  THEN
1088       IF ( myid == 0 )  THEN
1089          PRINT*, '+++ check_parameters: averaging_interval_pr=', &
1090                       averaging_interval_pr, ' must be <= dt_dopr=', dt_dopr
1091       ENDIF
1092       CALL local_stop
1093    ENDIF
1094
1095    IF ( averaging_interval_sp == 9999999.9 )  THEN
1096       averaging_interval_sp = averaging_interval
1097    ENDIF
1098
1099    IF ( averaging_interval_sp > dt_dosp )  THEN
1100       IF ( myid == 0 )  THEN
1101          PRINT*, '+++ check_parameters: averaging_interval_sp=', &
1102                       averaging_interval_sp, ' must be <= dt_dosp=', dt_dosp
1103       ENDIF
1104       CALL local_stop
1105    ENDIF
1106
1107!
1108!-- Set the default interval for profiles entering the temporal average
1109    IF ( dt_averaging_input_pr == 9999999.9 )  THEN
1110       dt_averaging_input_pr = dt_averaging_input
1111    ENDIF
1112
1113!
1114!-- Set the default interval for the output of timeseries to a reasonable
1115!-- value (tries to minimize the number of calls of flow_statistics)
1116    IF ( dt_dots == 9999999.9 )  THEN
1117       IF ( averaging_interval_pr == 0.0 )  THEN
1118          dt_dots = MIN( dt_run_control, dt_dopr )
1119       ELSE
1120          dt_dots = MIN( dt_run_control, dt_averaging_input_pr )
1121       ENDIF
1122    ENDIF
1123
1124!
1125!-- Check the sample rate for averaging (first for 3d-data, then for profiles)
1126    IF ( dt_averaging_input > averaging_interval )  THEN
1127       IF ( myid == 0 )  THEN
1128          PRINT*, '+++ check_parameters: dt_averaging_input=',                &
1129                       dt_averaging_input, ' must be <= averaging_interval=', &
1130                       averaging_interval
1131       ENDIF
1132       CALL local_stop
1133    ENDIF
1134
1135    IF ( dt_averaging_input_pr > averaging_interval_pr )  THEN
1136       IF ( myid == 0 )  THEN
1137          PRINT*, '+++ check_parameters: dt_averaging_input_pr=', &
1138                       dt_averaging_input_pr,                     &
1139                       ' must be <= averaging_interval_pr=',      &
1140                       averaging_interval_pr
1141       ENDIF
1142       CALL local_stop
1143    ENDIF
1144
1145!
1146!-- Determine the number of output profiles and check whether they are
1147!-- permissible
1148    DO  WHILE ( data_output_pr(dopr_n+1) /= '          ' )
1149
1150       dopr_n = dopr_n + 1
1151       i = dopr_n
1152
1153!
1154!--    Determine internal profile number (for hom, homs)
1155!--    and store height levels
1156       SELECT CASE ( TRIM( data_output_pr(i) ) )
1157
1158          CASE ( 'u', '#u' )
1159             dopr_index(i) = 1
1160             hom(:,2,1,:)  = SPREAD( zu, 2, statistic_regions+1 )
1161             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1162                dopr_initial_index(i) = 5
1163                hom(:,2,5,:)          = SPREAD( zu, 2, statistic_regions+1 )
1164                data_output_pr(i)     = data_output_pr(i)(2:)
1165             ENDIF
1166
1167          CASE ( 'v', '#v' )
1168             dopr_index(i) = 2
1169             hom(:,2,2,:) = SPREAD( zu, 2, statistic_regions+1 )
1170             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1171                dopr_initial_index(i) = 6
1172                hom(:,2,6,:)          = SPREAD( zu, 2, statistic_regions+1 )
1173                data_output_pr(i)     = data_output_pr(i)(2:)
1174             ENDIF
1175
1176          CASE ( 'w' )
1177             dopr_index(i) = 3
1178             hom(:,2,3,:) = SPREAD( zw, 2, statistic_regions+1 )
1179
1180          CASE ( 'pt', '#pt' )
1181             IF ( .NOT. cloud_physics ) THEN
1182                dopr_index(i) = 4
1183                hom(:,2,4,:)  = SPREAD( zu, 2, statistic_regions+1 )
1184                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1185                   dopr_initial_index(i) = 7
1186                   hom(:,2,7,:)          = SPREAD( zu, 2, statistic_regions+1 )
1187                   hom(nzb,2,7,:)        = 0.0    ! weil zu(nzb) negativ ist
1188                   data_output_pr(i)     = data_output_pr(i)(2:)
1189                ENDIF
1190             ELSE
1191                dopr_index(i) = 43
1192                hom(:,2,43,:)  = SPREAD( zu, 2, statistic_regions+1 )
1193                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1194                   dopr_initial_index(i) = 28
1195                   hom(:,2,28,:)         = SPREAD( zu, 2, statistic_regions+1 )
1196                   hom(nzb,2,28,:)       = 0.0    ! weil zu(nzb) negativ ist
1197                   data_output_pr(i)     = data_output_pr(i)(2:)
1198                ENDIF
1199             ENDIF
1200
1201          CASE ( 'e' )
1202             dopr_index(i)  = 8
1203             hom(:,2,8,:)   = SPREAD( zu, 2, statistic_regions+1 )
1204             hom(nzb,2,8,:) = 0.0
1205
1206          CASE ( 'km', '#km' )
1207             dopr_index(i)  = 9
1208             hom(:,2,9,:)   = SPREAD( zu, 2, statistic_regions+1 )
1209             hom(nzb,2,9,:) = 0.0
1210             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1211                dopr_initial_index(i) = 23
1212                hom(:,2,23,:)         = hom(:,2,9,:)
1213                data_output_pr(i)     = data_output_pr(i)(2:)
1214             ENDIF
1215
1216          CASE ( 'kh', '#kh' )
1217             dopr_index(i)   = 10
1218             hom(:,2,10,:)   = SPREAD( zu, 2, statistic_regions+1 )
1219             hom(nzb,2,10,:) = 0.0
1220             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1221                dopr_initial_index(i) = 24
1222                hom(:,2,24,:)         = hom(:,2,10,:)
1223                data_output_pr(i)     = data_output_pr(i)(2:)
1224             ENDIF
1225
1226          CASE ( 'l', '#l' )
1227             dopr_index(i)   = 11
1228             hom(:,2,11,:)   = SPREAD( zu, 2, statistic_regions+1 )
1229             hom(nzb,2,11,:) = 0.0
1230             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1231                dopr_initial_index(i) = 25
1232                hom(:,2,25,:)         = hom(:,2,11,:)
1233                data_output_pr(i)     = data_output_pr(i)(2:)
1234             ENDIF
1235
1236          CASE ( 'w"u"' )
1237             dopr_index(i) = 12
1238             hom(:,2,12,:) = SPREAD( zw, 2, statistic_regions+1 )
1239             IF ( prandtl_layer )  hom(nzb,2,12,:) = zu(1)
1240
1241          CASE ( 'w*u*' )
1242             dopr_index(i) = 13
1243             hom(:,2,13,:) = SPREAD( zw, 2, statistic_regions+1 )
1244
1245          CASE ( 'w"v"' )
1246             dopr_index(i) = 14
1247             hom(:,2,14,:) = SPREAD( zw, 2, statistic_regions+1 )
1248             IF ( prandtl_layer )  hom(nzb,2,14,:) = zu(1)
1249
1250          CASE ( 'w*v*' )
1251             dopr_index(i) = 15
1252             hom(:,2,15,:) = SPREAD( zw, 2, statistic_regions+1 )
1253
1254          CASE ( 'w"pt"' )
1255             dopr_index(i) = 16
1256             hom(:,2,16,:) = SPREAD( zw, 2, statistic_regions+1 )
1257
1258          CASE ( 'w*pt*' )
1259             dopr_index(i) = 17
1260             hom(:,2,17,:) = SPREAD( zw, 2, statistic_regions+1 )
1261
1262          CASE ( 'wpt' )
1263             dopr_index(i) = 18
1264             hom(:,2,18,:) = SPREAD( zw, 2, statistic_regions+1 )
1265
1266          CASE ( 'wu' )
1267             dopr_index(i) = 19
1268             hom(:,2,19,:) = SPREAD( zw, 2, statistic_regions+1 )
1269             IF ( prandtl_layer )  hom(nzb,2,19,:) = zu(1)
1270
1271          CASE ( 'wv' )
1272             dopr_index(i) = 20
1273             hom(:,2,20,:) = SPREAD( zw, 2, statistic_regions+1 )
1274             IF ( prandtl_layer )  hom(nzb,2,20,:) = zu(1)
1275
1276          CASE ( 'w*pt*BC' )
1277             dopr_index(i) = 21
1278             hom(:,2,21,:) = SPREAD( zw, 2, statistic_regions+1 )
1279
1280          CASE ( 'wptBC' )
1281             dopr_index(i) = 22
1282             hom(:,2,22,:) = SPREAD( zw, 2, statistic_regions+1 )
1283
1284          CASE ( 'u*2' )
1285             dopr_index(i) = 30
1286             hom(:,2,30,:) = SPREAD( zu, 2, statistic_regions+1 )
1287
1288          CASE ( 'v*2' )
1289             dopr_index(i) = 31
1290             hom(:,2,31,:) = SPREAD( zu, 2, statistic_regions+1 )
1291
1292          CASE ( 'w*2' )
1293             dopr_index(i) = 32
1294             hom(:,2,32,:) = SPREAD( zw, 2, statistic_regions+1 )
1295
1296          CASE ( 'pt*2' )
1297             dopr_index(i) = 33
1298             hom(:,2,33,:) = SPREAD( zu, 2, statistic_regions+1 )
1299
1300          CASE ( 'e*' )
1301             dopr_index(i) = 34
1302             hom(:,2,34,:) = SPREAD( zu, 2, statistic_regions+1 )
1303
1304          CASE ( 'w*2pt*' )
1305             dopr_index(i) = 35
1306             hom(:,2,35,:) = SPREAD( zw, 2, statistic_regions+1 )
1307
1308          CASE ( 'w*pt*2' )
1309             dopr_index(i) = 36
1310             hom(:,2,36,:) = SPREAD( zw, 2, statistic_regions+1 )
1311
1312          CASE ( 'w*e*' )
1313             dopr_index(i) = 37
1314             hom(:,2,37,:) = SPREAD( zw, 2, statistic_regions+1 )
1315
1316          CASE ( 'w*3' )
1317             dopr_index(i) = 38
1318             hom(:,2,38,:) = SPREAD( zw, 2, statistic_regions+1 )
1319
1320          CASE ( 'Sw' )
1321             dopr_index(i) = 39
1322             hom(:,2,39,:) = SPREAD( zw, 2, statistic_regions+1 )
1323
1324          CASE ( 'q', '#q' )
1325             IF ( .NOT. cloud_physics )  THEN
1326                IF ( myid == 0 )  THEN
1327                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1328                           data_output_pr(i),                          &
1329                           '    is not implemented for cloud_physics = FALSE'
1330                ENDIF
1331                CALL local_stop
1332             ELSE
1333                dopr_index(i) = 41
1334                hom(:,2,41,:)  = SPREAD( zu, 2, statistic_regions+1 )
1335                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1336                   dopr_initial_index(i) = 26
1337                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1338                   hom(nzb,2,26,:)       = 0.0    ! weil zu(nzb) negativ ist
1339                   data_output_pr(i)     = data_output_pr(i)(2:)
1340                ENDIF
1341             ENDIF
1342
1343          CASE ( 's', '#s' )
1344             IF ( .NOT. passive_scalar )  THEN
1345                IF ( myid == 0 )  THEN
1346                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1347                           data_output_pr(i),                          &
1348                           '    is not implemented for passive_scalar = FALSE'
1349                ENDIF
1350                CALL local_stop
1351             ELSE
1352                dopr_index(i) = 41
1353                hom(:,2,41,:)  = SPREAD( zu, 2, statistic_regions+1 )
1354                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1355                   dopr_initial_index(i) = 26
1356                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1357                   hom(nzb,2,26,:)       = 0.0    ! weil zu(nzb) negativ ist
1358                   data_output_pr(i)     = data_output_pr(i)(2:)
1359                ENDIF
1360             ENDIF
1361
1362          CASE ( 'qv', '#qv' )
1363             IF ( .NOT. cloud_physics ) THEN
1364                dopr_index(i) = 41
1365                hom(:,2,41,:)  = SPREAD( zu, 2, statistic_regions+1 )
1366                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1367                   dopr_initial_index(i) = 26
1368                   hom(:,2,26,:)         = SPREAD( zu, 2, statistic_regions+1 )
1369                   hom(nzb,2,26,:)       = 0.0    ! weil zu(nzb) negativ ist
1370                   data_output_pr(i)     = data_output_pr(i)(2:)
1371                ENDIF
1372             ELSE
1373                dopr_index(i) = 42
1374                hom(:,2,42,:)  = SPREAD( zu, 2, statistic_regions+1 )
1375                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1376                   dopr_initial_index(i) = 27
1377                   hom(:,2,27,:)         = SPREAD( zu, 2, statistic_regions+1 )
1378                   hom(nzb,2,27,:)       = 0.0    ! weil zu(nzb) negativ ist
1379                   data_output_pr(i)     = data_output_pr(i)(2:)
1380                ENDIF
1381             ENDIF
1382
1383          CASE ( 'lpt', '#lpt' )
1384             IF ( .NOT. cloud_physics ) THEN
1385                IF ( myid == 0 )  THEN
1386                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1387                           data_output_pr(i),                          &
1388                           '    is not implemented for cloud_physics = FALSE'
1389                ENDIF
1390                CALL local_stop
1391             ELSE
1392                dopr_index(i) = 4
1393                hom(:,2,4,:)  = SPREAD( zu, 2, statistic_regions+1 )
1394                IF ( data_output_pr(i)(1:1) == '#' )  THEN
1395                   dopr_initial_index(i) = 7
1396                   hom(:,2,7,:)          = SPREAD( zu, 2, statistic_regions+1 )
1397                   hom(nzb,2,7,:)        = 0.0    ! weil zu(nzb) negativ ist
1398                   data_output_pr(i)     = data_output_pr(i)(2:)
1399                ENDIF
1400             ENDIF
1401
1402          CASE ( 'vpt', '#vpt' )
1403             dopr_index(i) = 44
1404             hom(:,2,44,:)  = SPREAD( zu, 2, statistic_regions+1 )
1405             IF ( data_output_pr(i)(1:1) == '#' )  THEN
1406                dopr_initial_index(i) = 29
1407                hom(:,2,29,:)         = SPREAD( zu, 2, statistic_regions+1 )
1408                hom(nzb,2,29,:)       = 0.0    ! weil zu(nzb) negativ ist
1409                data_output_pr(i)     = data_output_pr(i)(2:)
1410             ENDIF
1411
1412          CASE ( 'w"vpt"' )
1413             dopr_index(i) = 45
1414             hom(:,2,45,:) = SPREAD( zw, 2, statistic_regions+1 )
1415
1416          CASE ( 'w*vpt*' )
1417             dopr_index(i) = 46
1418             hom(:,2,46,:) = SPREAD( zw, 2, statistic_regions+1 )
1419
1420          CASE ( 'wvpt' )
1421             dopr_index(i) = 47
1422             hom(:,2,47,:) = SPREAD( zw, 2, statistic_regions+1 )
1423
1424          CASE ( 'w"q"' )
1425             IF ( .NOT. cloud_physics ) THEN
1426                IF ( myid == 0 )  THEN
1427                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1428                           data_output_pr(i),                          &
1429                           '    is not implemented for cloud_physics = FALSE'
1430                ENDIF
1431                CALL local_stop
1432             ELSE
1433                dopr_index(i) = 48
1434                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
1435             ENDIF
1436
1437          CASE ( 'w*q*' )
1438             IF ( .NOT. cloud_physics ) THEN
1439                IF ( myid == 0 )  THEN
1440                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1441                           data_output_pr(i),                          &
1442                           '    is not implemented for cloud_physics = FALSE'
1443                ENDIF
1444                CALL local_stop
1445             ELSE
1446                dopr_index(i) = 49
1447                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
1448             ENDIF
1449
1450          CASE ( 'wq' )
1451             IF ( .NOT. cloud_physics ) THEN
1452                IF ( myid == 0 )  THEN
1453                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1454                           data_output_pr(i),                          &
1455                           '    is not implemented for cloud_physics = FALSE'
1456                ENDIF
1457                CALL local_stop
1458             ELSE
1459                dopr_index(i) = 50
1460                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
1461             ENDIF
1462
1463          CASE ( 'w"s"' )
1464             IF ( .NOT. passive_scalar ) THEN
1465                IF ( myid == 0 )  THEN
1466                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1467                           data_output_pr(i),                          &
1468                           '    is not implemented for passive_scalar = FALSE'
1469                ENDIF
1470                CALL local_stop
1471             ELSE
1472                dopr_index(i) = 48
1473                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
1474             ENDIF
1475
1476          CASE ( 'w*s*' )
1477             IF ( .NOT. passive_scalar ) THEN
1478                IF ( myid == 0 )  THEN
1479                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1480                           data_output_pr(i),                          &
1481                           '    is not implemented for passive_scalar = FALSE'
1482                ENDIF
1483                CALL local_stop
1484             ELSE
1485                dopr_index(i) = 49
1486                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
1487             ENDIF
1488
1489          CASE ( 'ws' )
1490             IF ( .NOT. passive_scalar ) THEN
1491                IF ( myid == 0 )  THEN
1492                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1493                           data_output_pr(i),                          &
1494                           '    is not implemented for passive_scalar = FALSE'
1495                ENDIF
1496                CALL local_stop
1497             ELSE
1498                dopr_index(i) = 50
1499                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
1500             ENDIF
1501
1502          CASE ( 'w"qv"' )
1503             IF ( moisture  .AND.  .NOT. cloud_physics ) &
1504             THEN
1505                dopr_index(i) = 48
1506                hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 )
1507             ELSEIF( moisture .AND. cloud_physics ) THEN
1508                dopr_index(i) = 51
1509                hom(:,2,51,:) = SPREAD( zw, 2, statistic_regions+1 )
1510             ELSE
1511                IF ( myid == 0 )  THEN
1512                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1513                           data_output_pr(i),                          &
1514                           '    is not implemented for cloud_physics = FALSE', &
1515                           '    and                    moisture      = FALSE'
1516                ENDIF
1517                CALL local_stop                   
1518             ENDIF
1519
1520          CASE ( 'w*qv*' )
1521             IF ( moisture  .AND.  .NOT. cloud_physics ) &
1522             THEN
1523                dopr_index(i) = 49
1524                hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 )
1525             ELSEIF( moisture .AND. cloud_physics ) THEN
1526                dopr_index(i) = 52
1527                hom(:,2,52,:) = SPREAD( zw, 2, statistic_regions+1 )
1528             ELSE
1529                IF ( myid == 0 )  THEN
1530                   PRINT*, '+++ check_parameters:  data_output_pr = ',         &
1531                           data_output_pr(i),                                  &
1532                           '    is not implemented for cloud_physics = FALSE', &
1533                           '                       and moisture      = FALSE'
1534                ENDIF
1535                CALL local_stop                   
1536             ENDIF
1537
1538          CASE ( 'wqv' )
1539             IF ( moisture  .AND.  .NOT. cloud_physics ) &
1540             THEN
1541                dopr_index(i) = 50
1542                hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 )
1543             ELSEIF( moisture .AND. cloud_physics ) THEN
1544                dopr_index(i) = 53
1545                hom(:,2,53,:) = SPREAD( zw, 2, statistic_regions+1 )
1546             ELSE
1547                IF ( myid == 0 )  THEN
1548                   PRINT*, '+++ check_parameters:  data_output_pr = ',         &
1549                           data_output_pr(i),                                  &
1550                           '    is not implemented for cloud_physics = FALSE', &
1551                           '                       and moisture      = FALSE'
1552                ENDIF
1553                CALL local_stop                   
1554             ENDIF
1555
1556          CASE ( 'ql' )
1557             IF ( .NOT. cloud_physics  .AND.  .NOT. cloud_droplets )  THEN
1558                IF ( myid == 0 )  THEN
1559                   PRINT*, '+++ check_parameters:  data_output_pr = ', &
1560                           data_output_pr(i),                          &
1561                           '    is not implemented for cloud_physics = FALSE'
1562                ENDIF
1563                CALL local_stop
1564             ELSE
1565                dopr_index(i) = 54
1566                hom(:,2,54,:)  = SPREAD( zu, 2, statistic_regions+1 )
1567             ENDIF
1568
1569          CASE ( 'w*u*u*/dz' )
1570             dopr_index(i) = 55
1571             hom(:,2,55,:) = SPREAD( zu, 2, statistic_regions+1 )
1572
1573          CASE ( 'w*p*/dz' )
1574             dopr_index(i) = 56
1575             hom(:,2,56,:) = SPREAD( zu, 2, statistic_regions+1 )
1576
1577          CASE ( 'w"e/dz' )
1578             dopr_index(i) = 57
1579             hom(:,2,57,:) = SPREAD( zu, 2, statistic_regions+1 )
1580
1581          CASE ( 'u"pt"' )
1582             dopr_index(i) = 58
1583             hom(:,2,58,:) = SPREAD( zu, 2, statistic_regions+1 )
1584
1585          CASE ( 'u*pt*' )
1586             dopr_index(i) = 59
1587             hom(:,2,59,:) = SPREAD( zu, 2, statistic_regions+1 )
1588
1589          CASE ( 'upt_t' )
1590             dopr_index(i) = 60
1591             hom(:,2,60,:) = SPREAD( zu, 2, statistic_regions+1 )
1592
1593          CASE ( 'v"pt"' )
1594             dopr_index(i) = 61
1595             hom(:,2,61,:) = SPREAD( zu, 2, statistic_regions+1 )
1596             
1597          CASE ( 'v*pt*' )
1598             dopr_index(i) = 62
1599             hom(:,2,62,:) = SPREAD( zu, 2, statistic_regions+1 )
1600
1601          CASE ( 'vpt_t' )
1602             dopr_index(i) = 63
1603             hom(:,2,63,:) = SPREAD( zu, 2, statistic_regions+1 )
1604
1605
1606          CASE DEFAULT
1607             IF ( myid == 0 )  THEN
1608                PRINT*, '+++ check_parameters:  unknown output profile:  ', &
1609                        'data_output_pr = ', data_output_pr(i)
1610             ENDIF
1611             CALL local_stop
1612
1613       END SELECT
1614!
1615!--    Check to which of the predefined coordinate systems the profile belongs
1616       DO  k = 1, crmax
1617          IF ( INDEX( cross_profiles(k), ' '//TRIM( data_output_pr(i) )//' ' ) &
1618               /=0 ) &
1619          THEN
1620             dopr_crossindex(i) = k
1621             EXIT
1622          ENDIF
1623       ENDDO
1624!
1625!--    Generate the text for the labels of the PROFIL output file. "-characters
1626!--    must be substituted, otherwise PROFIL would interpret them as TeX
1627!--    control characters
1628       dopr_label(i) = data_output_pr(i)
1629       position = INDEX( dopr_label(i) , '"' )
1630       DO WHILE ( position /= 0 )
1631          dopr_label(i)(position:position) = ''''
1632          position = INDEX( dopr_label(i) , '"' )
1633       ENDDO
1634
1635    ENDDO
1636
1637!
1638!-- y-value range of the coordinate system (PROFIL).
1639!-- x-value range determined in plot_1d.
1640    cross_uymin = 0.0
1641    IF ( z_max_do1d == -1.0 )  THEN
1642       cross_uymax = zu(nzt+1)
1643    ELSEIF ( z_max_do1d < zu(nzb+1)  .OR.  z_max_do1d > zu(nzt+1) )  THEN
1644       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  z_max_do1d=',  &
1645                                 z_max_do1d,' must be >= ', zu(nzb+1),  &
1646                                 ' or <= ', zu(nzt+1)
1647       CALL local_stop
1648    ELSE
1649       cross_uymax = z_max_do1d
1650    ENDIF
1651
1652!
1653!-- Check whether the chosen normalizing factor for the coordinate systems is
1654!-- permissible
1655    DO  i = 1, crmax
1656       SELECT CASE ( TRIM( cross_normalized_x(i) ) )  ! TRIM required on IBM
1657
1658          CASE ( '', 'wpt0', 'ws2', 'tsw2', 'ws3', 'ws2tsw', 'wstsw2' )
1659             j = 0
1660
1661          CASE DEFAULT
1662             IF ( myid == 0 )  THEN
1663                PRINT*, '+++ check_parameters: unknown normalize method'
1664                PRINT*, '    cross_normalized_x="',cross_normalized_x(i),'"'
1665             ENDIF
1666             CALL local_stop
1667
1668       END SELECT
1669       SELECT CASE ( TRIM( cross_normalized_y(i) ) )  ! TRIM required on IBM
1670
1671          CASE ( '', 'z_i' )
1672             j = 0
1673
1674          CASE DEFAULT
1675             IF ( myid == 0 )  THEN
1676                PRINT*, '+++ check_parameters: unknown normalize method'
1677                PRINT*, '    cross_normalized_y="',cross_normalized_y(i),'"'
1678             ENDIF
1679             CALL local_stop
1680
1681       END SELECT
1682    ENDDO
1683!
1684!-- Check normalized y-value range of the coordinate system (PROFIL)
1685    IF ( z_max_do1d_normalized /= -1.0  .AND.  z_max_do1d_normalized <= 0.0 ) &
1686    THEN
1687       IF ( myid == 0 )  PRINT*,'+++ check_parameters:  z_max_do1d_normalize', &
1688                                'd=', z_max_do1d_normalized, ' must be >= 0.0'
1689       CALL local_stop
1690    ENDIF
1691
1692
1693!
1694!-- Append user-defined data output variables to the standard data output
1695    IF ( data_output_user(1) /= ' ' )  THEN
1696       i = 1
1697       DO  WHILE ( data_output(i) /= ' '  .AND.  i <= 100 )
1698          i = i + 1
1699       ENDDO
1700       j = 1
1701       DO  WHILE ( data_output_user(j) /= ' '  .AND.  j <= 100 )
1702          IF ( i > 100 )  THEN
1703             IF ( myid == 0 )  THEN
1704                PRINT*, '+++ check_parameters: number of output quantitities', &
1705                             ' given by data_output and data_output_user'
1706                PRINT*, '                      exceeds the limit of 100'
1707             ENDIF
1708             CALL local_stop
1709          ENDIF
1710          data_output(i) = data_output_user(j)
1711          i = i + 1
1712          j = j + 1
1713       ENDDO
1714    ENDIF
1715
1716!
1717!-- Check and set steering parameters for 2d/3d data output and averaging
1718    i   = 1
1719    DO  WHILE ( data_output(i) /= ' '  .AND.  i <= 100 )
1720!
1721!--    Check for data averaging
1722       ilen = LEN_TRIM( data_output(i) )
1723       j = 0                                                 ! no data averaging
1724       IF ( ilen > 3 )  THEN
1725          IF ( data_output(i)(ilen-2:ilen) == '_av' )  THEN
1726             j = 1                                           ! data averaging
1727             data_output(i) = data_output(i)(1:ilen-3)
1728          ENDIF
1729       ENDIF
1730!
1731!--    Check for cross section or volume data
1732       ilen = LEN_TRIM( data_output(i) )
1733       k = 0                                                   ! 3d data
1734       var = data_output(i)(1:ilen)
1735       IF ( ilen > 3 )  THEN
1736          IF ( data_output(i)(ilen-2:ilen) == '_xy'  .OR. &
1737               data_output(i)(ilen-2:ilen) == '_xz'  .OR. &
1738               data_output(i)(ilen-2:ilen) == '_yz' )  THEN
1739             k = 1                                             ! 2d data
1740             var = data_output(i)(1:ilen-3)
1741          ENDIF
1742       ENDIF
1743!
1744!--    Check for allowed value and set units
1745       SELECT CASE ( TRIM( var ) )
1746
1747          CASE ( 'e' )
1748             IF ( constant_diffusion )  THEN
1749                IF ( myid == 0 )  THEN
1750                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1751                                '" requires constant_diffusion = .FALSE.'
1752                ENDIF
1753                CALL local_stop
1754             ENDIF
1755             unit = 'm2/s2'
1756
1757          CASE ( 'pc', 'pr' )
1758             IF ( .NOT. particle_advection )  THEN
1759                IF ( myid == 0 )  THEN
1760                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1761                                '" requires particle package'
1762                   PRINT*, '                      (mrun-option "-p particles")'
1763                ENDIF
1764                CALL local_stop
1765             ENDIF
1766             IF ( TRIM( var ) == 'pc' )  unit = 'number'
1767             IF ( TRIM( var ) == 'pr' )  unit = 'm'
1768
1769          CASE ( 'q', 'vpt' )
1770             IF ( .NOT. moisture )  THEN
1771                IF ( myid == 0 )  THEN
1772                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1773                                '" requires moisture = .TRUE.'
1774                ENDIF
1775                CALL local_stop
1776             ENDIF
1777             IF ( TRIM( var ) == 'q'   )  unit = 'kg/kg'
1778             IF ( TRIM( var ) == 'vpt' )  unit = 'K'
1779
1780          CASE ( 'ql' )
1781             IF ( .NOT. ( cloud_physics  .OR.  cloud_droplets ) )  THEN
1782                IF ( myid == 0 )  THEN
1783                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1784                                '" requires cloud_physics = .TRUE.'
1785                   PRINT*, '                      or cloud_droplets = .TRUE.'
1786                ENDIF
1787                CALL local_stop
1788             ENDIF
1789             unit = 'kg/kg'
1790
1791          CASE ( 'ql_c', 'ql_v', 'ql_vp' )
1792             IF ( .NOT. cloud_droplets )  THEN
1793                IF ( myid == 0 )  THEN
1794                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1795                                '" requires cloud_droplets = .TRUE.'
1796                ENDIF
1797                CALL local_stop
1798             ENDIF
1799             IF ( TRIM( var ) == 'ql_c'  )  unit = 'kg/kg'
1800             IF ( TRIM( var ) == 'ql_v'  )  unit = 'm3'
1801             IF ( TRIM( var ) == 'ql_vp' )  unit = 'none'
1802
1803          CASE ( 'qv' )
1804             IF ( .NOT. cloud_physics )  THEN
1805                IF ( myid == 0 )  THEN
1806                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1807                                '" requires cloud_physics = .TRUE.'
1808                ENDIF
1809                CALL local_stop
1810             ENDIF
1811             unit = 'kg/kg'
1812
1813          CASE ( 's' )
1814             IF ( .NOT. passive_scalar )  THEN
1815                IF ( myid == 0 )  THEN
1816                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1817                                '" requires passive_scalar = .TRUE.'
1818                ENDIF
1819                CALL local_stop
1820             ENDIF
1821             unit = 'conc'
1822
1823          CASE ( 'u*', 't*', 'lwp*' )
1824             IF ( k == 0  .OR.  data_output(i)(ilen-2:ilen) /= '_xy' )  THEN
1825                IF ( myid == 0 )  THEN
1826                   PRINT*, '+++ check_parameters:  illegal value for data_',&
1827                                'output: "', TRIM( var ), '" is only allowed'
1828                   PRINT*, '                       for horizontal cross section'
1829                ENDIF
1830                CALL local_stop
1831             ENDIF
1832             IF ( TRIM( var ) == 'lwp*'  .AND.  .NOT. cloud_physics )  THEN
1833                IF ( myid == 0 )  THEN
1834                   PRINT*, '+++ check_parameters: output of "', TRIM( var ), &
1835                                '" requires cloud_physics = .TRUE.'
1836                ENDIF
1837                CALL local_stop
1838             ENDIF
1839             IF ( TRIM( var ) == 'u*'   )  unit = 'm/s'
1840             IF ( TRIM( var ) == 't*'   )  unit = 'K'
1841             IF ( TRIM( var ) == 'lwp*' )  unit = 'kg/kg*m'
1842
1843          CASE ( 'p', 'pt', 'u', 'v', 'w' )
1844             IF ( TRIM( var ) == 'p'  )  unit = 'Pa'
1845             IF ( TRIM( var ) == 'pt' )  unit = 'K'
1846             IF ( TRIM( var ) == 'u'  )  unit = 'm/s'
1847             IF ( TRIM( var ) == 'v'  )  unit = 'm/s'
1848             IF ( TRIM( var ) == 'w'  )  unit = 'm/s'
1849             CONTINUE
1850
1851          CASE DEFAULT
1852             CALL user_check_data_output( var, unit )
1853
1854             IF ( unit == 'illegal' )  THEN
1855                IF ( myid == 0 )  THEN
1856                   IF ( data_output_user(1) /= ' ' )  THEN
1857                      PRINT*, '+++ check_parameters:  illegal value for data_',&
1858                                   'output or data_output_user: "',            &
1859                                   TRIM( data_output(i) ), '"'
1860                   ELSE
1861                      PRINT*, '+++ check_parameters:  illegal value for data_',&
1862                                   'output: "', TRIM( data_output(i) ), '"'
1863                   ENDIF
1864                ENDIF
1865                CALL local_stop
1866             ENDIF
1867
1868       END SELECT
1869!
1870!--    Set the internal steering parameters appropriately
1871       IF ( k == 0 )  THEN
1872          do3d_no(j)              = do3d_no(j) + 1
1873          do3d(j,do3d_no(j))      = data_output(i)
1874          do3d_unit(j,do3d_no(j)) = unit
1875       ELSE
1876          do2d_no(j)              = do2d_no(j) + 1
1877          do2d(j,do2d_no(j))      = data_output(i)
1878          do2d_unit(j,do2d_no(j)) = unit
1879          IF ( data_output(i)(ilen-2:ilen) == '_xy' )  THEN
1880             data_output_xy(j) = .TRUE.
1881          ENDIF
1882          IF ( data_output(i)(ilen-2:ilen) == '_xz' )  THEN
1883             data_output_xz(j) = .TRUE.
1884          ENDIF
1885          IF ( data_output(i)(ilen-2:ilen) == '_yz' )  THEN
1886             data_output_yz(j) = .TRUE.
1887          ENDIF
1888       ENDIF
1889
1890       IF ( j == 1 )  THEN
1891!
1892!--       Check, if variable is already subject to averaging
1893          found = .FALSE.
1894          DO  k = 1, doav_n
1895             IF ( TRIM( doav(k) ) == TRIM( var ) )  found = .TRUE.
1896          ENDDO
1897
1898          IF ( .NOT. found )  THEN
1899             doav_n = doav_n + 1
1900             doav(doav_n) = var
1901          ENDIF
1902       ENDIF
1903
1904       i = i + 1
1905    ENDDO
1906
1907!
1908!-- Store sectional planes in one shared array
1909    section(:,1) = section_xy
1910    section(:,2) = section_xz
1911    section(:,3) = section_yz
1912
1913!
1914!-- Upper plot limit (grid point value) for 1D profiles
1915    IF ( z_max_do1d == -1.0 )  THEN
1916       nz_do1d = nzt+1
1917    ELSE
1918       DO  k = nzb+1, nzt+1
1919          nz_do1d = k
1920          IF ( zw(k) > z_max_do1d )  EXIT
1921       ENDDO
1922    ENDIF
1923
1924!
1925!-- Upper plot limit for 2D vertical sections
1926    IF ( z_max_do2d == -1.0 )  z_max_do2d = zu(nzt)
1927    IF ( z_max_do2d < zu(nzb+1)  .OR.  z_max_do2d > zu(nzt) )  THEN
1928       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  z_max_do2d=',        &
1929                                 z_max_do2d, ' must be >= ', zu(nzb+1),       &
1930                                 '(zu(nzb+1)) and <= ', zu(nzt), ' (zu(nzt))'
1931       CALL local_stop
1932    ENDIF
1933
1934!
1935!-- Upper plot limit for 3D arrays
1936    IF ( nz_do3d == -9999 )  nz_do3d = nzt + 1
1937
1938!
1939!-- Determine and check accuracy for compressed 3D plot output
1940    IF ( do3d_compress )  THEN
1941!
1942!--    Compression only permissible on T3E machines
1943       IF ( host(1:3) /= 't3e' )  THEN
1944          IF ( myid == 0 )  THEN
1945             PRINT*, '+++ check_parameters: do3d_compress = .TRUE. not allow', &
1946                          'ed on host "', TRIM( host ), '"'
1947          ENDIF
1948          CALL local_stop
1949       ENDIF
1950
1951       i = 1
1952       DO  WHILE ( do3d_comp_prec(i) /= ' ' )
1953
1954          ilen = LEN_TRIM( do3d_comp_prec(i) )
1955          IF ( LLT( do3d_comp_prec(i)(ilen:ilen), '0' ) .OR. &
1956               LGT( do3d_comp_prec(i)(ilen:ilen), '9' ) )  THEN
1957             IF ( myid == 0 )  THEN
1958                PRINT*, '+++ check_parameters: illegal precision: ', &
1959                        'do3d_comp_prec(', i, ')="', TRIM(do3d_comp_prec(i)),'"'
1960             ENDIF
1961             CALL local_stop
1962          ENDIF
1963
1964          prec = IACHAR( do3d_comp_prec(i)(ilen:ilen) ) - IACHAR( '0' )
1965          var = do3d_comp_prec(i)(1:ilen-1)
1966
1967          SELECT CASE ( var )
1968
1969             CASE ( 'u' )
1970                j = 1
1971             CASE ( 'v' )
1972                j = 2
1973             CASE ( 'w' )
1974                j = 3
1975             CASE ( 'p' )
1976                j = 4
1977             CASE ( 'pt' )
1978                j = 5
1979
1980             CASE DEFAULT
1981                IF ( myid == 0 )  THEN
1982                   PRINT*, '+++ check_parameters: unknown variable in ', &
1983                               'assignment'
1984                   PRINT*, '    do3d_comp_prec(', i, ')="', &
1985                           TRIM( do3d_comp_prec(i) ),'"'
1986                ENDIF
1987                CALL local_stop               
1988
1989          END SELECT
1990
1991          plot_3d_precision(j)%precision = prec
1992          i = i + 1
1993
1994       ENDDO
1995    ENDIF
1996
1997!
1998!-- Check the data output format(s)
1999    IF ( data_output_format(1) == ' ' )  THEN
2000!
2001!--    Default value
2002       netcdf_output = .TRUE.
2003    ELSE
2004       i = 1
2005       DO  WHILE ( data_output_format(i) /= ' ' )
2006
2007          SELECT CASE ( data_output_format(i) )
2008
2009             CASE ( 'netcdf' )
2010                netcdf_output = .TRUE.
2011             CASE ( 'iso2d' )
2012                iso2d_output  = .TRUE.
2013             CASE ( 'profil' )
2014                profil_output = .TRUE.
2015             CASE ( 'avs' )
2016                avs_output    = .TRUE.
2017
2018             CASE DEFAULT
2019                IF ( myid == 0 )  THEN
2020                   PRINT*, '+++ check_parameters:'
2021                   PRINT*, '    unknown value for data_output_format "', &
2022                                TRIM( data_output_format(i) ),'"'
2023                ENDIF
2024                CALL local_stop               
2025
2026          END SELECT
2027
2028          i = i + 1
2029          IF ( i > 10 )  EXIT
2030
2031       ENDDO
2032
2033    ENDIF
2034
2035!
2036!-- Check netcdf precison
2037    ldum = .FALSE.
2038    CALL define_netcdf_header( 'ch', ldum, 0 )
2039
2040!
2041!-- Check, whether a constant diffusion coefficient shall be used
2042    IF ( km_constant /= -1.0 )  THEN
2043       IF ( km_constant < 0.0 )  THEN
2044          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  km_constant=', &
2045                                    km_constant, ' < 0.0'
2046          CALL local_stop
2047       ELSE
2048          IF ( prandtl_number < 0.0 )  THEN
2049             IF ( myid == 0 )  PRINT*,'+++ check_parameters:  prandtl_number=',&
2050                                      prandtl_number, ' < 0.0'
2051             CALL local_stop
2052          ENDIF
2053          constant_diffusion = .TRUE.
2054
2055          IF ( prandtl_layer )  THEN
2056             IF ( myid == 0 )  PRINT*, '+++ check_parameters:  prandtl_layer ',&
2057                                       'is not allowed with fixed value of km'
2058             CALL local_stop
2059          ENDIF
2060       ENDIF
2061    ENDIF
2062
2063!
2064!-- In case of non-cyclic lateral boundaries, set the default maximum value
2065!-- for the horizontal diffusivity used within the outflow damping layer,
2066!-- and check/set the width of the damping layer
2067    IF ( bc_lr /= 'cyclic' ) THEN
2068       IF ( km_damp_max == -1.0 )  THEN
2069          km_damp_max = 0.5 * dx
2070       ENDIF
2071       IF ( outflow_damping_width == -1.0 )  THEN
2072          outflow_damping_width = MIN( 20, nx/2 )
2073       ENDIF
2074       IF ( outflow_damping_width <= 0  .OR.  outflow_damping_width > nx )  THEN
2075          IF ( myid == 0 )  PRINT*, '+++ check_parameters: outflow_damping w',&
2076                                    'idth out of range'
2077          CALL local_stop
2078       ENDIF
2079    ENDIF
2080
2081    IF ( bc_ns /= 'cyclic' )  THEN
2082       IF ( km_damp_max == -1.0 )  THEN
2083          km_damp_max = 0.5 * dy
2084       ENDIF
2085       IF ( outflow_damping_width == -1.0 )  THEN
2086          outflow_damping_width = MIN( 20, ny/2 )
2087       ENDIF
2088       IF ( outflow_damping_width <= 0  .OR.  outflow_damping_width > ny )  THEN
2089          IF ( myid == 0 )  PRINT*, '+++ check_parameters: outflow_damping w',&
2090                                    'idth out of range'
2091          CALL local_stop
2092       ENDIF
2093    ENDIF
2094
2095!
2096!-- Check value range for rif
2097    IF ( rif_min >= rif_max )  THEN
2098       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  rif_min=', rif_min, &
2099                                 ' must be less than rif_max=', rif_max
2100       CALL local_stop
2101    ENDIF
2102
2103!
2104!-- Determine upper and lower hight level indices for random perturbations
2105    IF ( disturbance_level_b == -1.0 )  THEN
2106       disturbance_level_b     = zu(nzb+3)
2107       disturbance_level_ind_b = nzb + 3
2108    ELSEIF ( disturbance_level_b < zu(3) )  THEN
2109       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_b=',&
2110                                 disturbance_level_b, ' must be >= ',zu(3),    &
2111                                 '(zu(3))'
2112       CALL local_stop
2113    ELSEIF ( disturbance_level_b > zu(nzt-2) )  THEN
2114       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_b=',&
2115                                 disturbance_level_b, ' must be <= ',zu(nzt-2),&
2116                                 '(zu(nzt-2))'
2117       CALL local_stop
2118    ELSE
2119       DO  k = 3, nzt-2
2120          IF ( disturbance_level_b <= zu(k) )  THEN
2121             disturbance_level_ind_b = k
2122             EXIT
2123          ENDIF
2124       ENDDO
2125    ENDIF
2126
2127    IF ( disturbance_level_t == -1.0 )  THEN
2128       disturbance_level_t     = zu(nzt/3)
2129       disturbance_level_ind_t = nzt / 3
2130    ELSEIF ( disturbance_level_t > zu(nzt-2) )  THEN
2131       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_t=',&
2132                                 disturbance_level_t, ' must be <= ',zu(nzt-2),&
2133                                 '(zu(nzt-2))'
2134       CALL local_stop
2135    ELSEIF ( disturbance_level_t < disturbance_level_b )  THEN
2136       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  disturbance_level_t=',&
2137                                 disturbance_level_t, ' must be >= ',          &
2138                                 'disturbance_level_b=', disturbance_level_b
2139       CALL local_stop
2140    ELSE
2141       DO  k = 3, nzt-2
2142          IF ( disturbance_level_t <= zu(k) )  THEN
2143             disturbance_level_ind_t = k
2144             EXIT
2145          ENDIF
2146       ENDDO
2147    ENDIF
2148
2149!
2150!-- Check again whether the levels determined this way are ok.
2151!-- Error may occur at automatic determination and too few grid points in
2152!-- z-direction.
2153    IF ( disturbance_level_ind_t < disturbance_level_ind_b )  THEN
2154       IF ( myid == 0 )  PRINT*, '+++ check_parameters:  ',                &
2155                                 'disturbance_level_ind_t=',               &
2156                                 disturbance_level_ind_t, ' must be >= ',  &
2157                                 'disturbance_level_ind_b=',               &
2158                                 disturbance_level_ind_b
2159       CALL local_stop
2160    ENDIF
2161
2162!
2163!-- Determine the horizontal index range for random perturbations.
2164!-- In case of non-cyclic horizontal boundaries, no perturbations are imposed
2165!-- near the inflow and the perturbation area is further limited to ...(1)
2166!-- after the initial phase of the flow.
2167    dist_nxl = 0;  dist_nxr = nx
2168    dist_nys = 0;  dist_nyn = ny
2169    IF ( bc_lr /= 'cyclic' )  THEN
2170       IF ( inflow_disturbance_begin == -1 )  THEN
2171          inflow_disturbance_begin = MIN( 10, nx/2 )
2172       ENDIF
2173       IF ( inflow_disturbance_begin < 0  .OR.  inflow_disturbance_begin > nx )&
2174       THEN
2175          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2176                                    '_begin out of range'
2177          CALL local_stop
2178       ENDIF
2179       IF ( inflow_disturbance_end == -1 )  THEN
2180          inflow_disturbance_end = MIN( 100, 3*nx/4 )
2181       ENDIF
2182       IF ( inflow_disturbance_end < 0  .OR.  inflow_disturbance_end > nx )    &
2183       THEN
2184          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2185                                    '_end out of range'
2186          CALL local_stop
2187       ENDIF
2188    ELSEIF ( bc_ns /= 'cyclic' )  THEN
2189       IF ( inflow_disturbance_begin == -1 )  THEN
2190          inflow_disturbance_begin = MIN( 10, ny/2 )
2191       ENDIF
2192       IF ( inflow_disturbance_begin < 0  .OR.  inflow_disturbance_begin > ny )&
2193       THEN
2194          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2195                                    '_begin out of range'
2196          CALL local_stop
2197       ENDIF
2198       IF ( inflow_disturbance_end == -1 )  THEN
2199          inflow_disturbance_end = MIN( 100, 3*ny/4 )
2200       ENDIF
2201       IF ( inflow_disturbance_end < 0  .OR.  inflow_disturbance_end > ny )    &
2202       THEN
2203          IF ( myid == 0 )  PRINT*, '+++ check_parameters: inflow_disturbance',&
2204                                    '_end out of range'
2205          CALL local_stop
2206       ENDIF
2207    ENDIF
2208
2209    IF ( outflow_l )  THEN
2210       dist_nxr    = nx - inflow_disturbance_begin
2211       dist_nxl(1) = nx - inflow_disturbance_end
2212    ELSEIF ( outflow_r )  THEN
2213       dist_nxl    = inflow_disturbance_begin
2214       dist_nxr(1) = inflow_disturbance_end
2215    ELSEIF ( outflow_s )  THEN
2216       dist_nyn    = ny - inflow_disturbance_begin
2217       dist_nys(1) = ny - inflow_disturbance_end
2218    ELSEIF ( outflow_n )  THEN
2219       dist_nys    = inflow_disturbance_begin
2220       dist_nyn(1) = inflow_disturbance_end
2221    ENDIF
2222
2223!
2224!-- Check random generator
2225    IF ( random_generator /= 'system-specific'  .AND. &
2226         random_generator /= 'numerical-recipes' )  THEN
2227       IF ( myid == 0 )  THEN
2228          PRINT*, '+++ check_parameters:'
2229          PRINT*, '    unknown random generator: random_generator=', &
2230                  random_generator
2231       ENDIF
2232       CALL local_stop
2233    ENDIF
2234
2235!
2236!-- Determine damping level index for 1D model
2237    IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
2238       IF ( damp_level_1d == -1.0 )  THEN
2239          damp_level_1d     = zu(nzt+1)
2240          damp_level_ind_1d = nzt + 1
2241       ELSEIF ( damp_level_1d < 0.0  .OR.  damp_level_1d > zu(nzt+1) )  THEN
2242          IF ( myid == 0 )  PRINT*, '+++ check_parameters:  damp_level_1d=', &
2243                                    damp_level_1d, ' must be > 0.0 and < ',  &
2244                                    'zu(nzt+1)', zu(nzt+1)
2245          CALL local_stop
2246       ELSE
2247          DO  k = 1, nzt+1
2248             IF ( damp_level_1d <= zu(k) )  THEN
2249                damp_level_ind_1d = k
2250                EXIT
2251             ENDIF
2252          ENDDO
2253       ENDIF
2254    ENDIF
2255!
2256!-- Check some other 1d-model parameters
2257    IF ( TRIM( mixing_length_1d ) /= 'as_in_3d_model'  .AND. &
2258         TRIM( mixing_length_1d ) /= 'blackadar' )  THEN
2259       IF ( myid == 0 )  PRINT*, '+++ check_parameters: mixing_length_1d = "', &
2260                                 TRIM( mixing_length_1d ), '" is unknown'
2261       CALL local_stop
2262    ENDIF
2263    IF ( TRIM( dissipation_1d ) /= 'as_in_3d_model'  .AND. &
2264         TRIM( dissipation_1d ) /= 'detering' )  THEN
2265       IF ( myid == 0 )  PRINT*, '+++ check_parameters: dissipation_1d = "', &
2266                                 TRIM( dissipation_1d ), '" is unknown'
2267       CALL local_stop
2268    ENDIF
2269
2270!
2271!-- Set time for the next user defined restart (time_restart is the
2272!-- internal parameter for steering restart events)
2273    IF ( restart_time /= 9999999.9 )  THEN
2274       IF ( restart_time > simulated_time )  time_restart = restart_time
2275    ELSE
2276!
2277!--    In case of a restart run, set internal parameter to default (no restart)
2278!--    if the NAMELIST-parameter restart_time is omitted
2279       time_restart = 9999999.9
2280    ENDIF
2281
2282!
2283!-- Set default value of the time needed to terminate a model run
2284    IF ( termination_time_needed == -1.0 )  THEN
2285       IF ( host(1:3) == 'ibm' )  THEN
2286          termination_time_needed = 300.0
2287       ELSE
2288          termination_time_needed = 35.0
2289       ENDIF
2290    ENDIF
2291
2292!
2293!-- Check the time needed to terminate a model run
2294    IF ( host(1:3) == 't3e' )  THEN
2295!
2296!--    Time needed must be at least 30 seconds on all CRAY machines, because
2297!--    MPP_TREMAIN gives the remaining CPU time only in steps of 30 seconds
2298       IF ( termination_time_needed <= 30.0 )  THEN
2299          IF ( myid == 0 )  THEN
2300             PRINT*, '+++ check_parameters:  termination_time_needed', &
2301                      termination_time_needed
2302             PRINT*, '                       must be > 30.0 on host "', host, &
2303                     '"'
2304          ENDIF
2305          CALL local_stop
2306       ENDIF
2307    ELSEIF ( host(1:3) == 'ibm' )  THEN
2308!
2309!--    On IBM-regatta machines the time should be at least 300 seconds,
2310!--    because the job time consumed before executing palm (for compiling,
2311!--    copying of files, etc.) has to be regarded
2312       IF ( termination_time_needed < 300.0 )  THEN
2313          IF ( myid == 0 )  THEN
2314             PRINT*, '+++ WARNING: check_parameters:  termination_time_',  &
2315                         'needed', termination_time_needed
2316             PRINT*, '                                should be >= 300.0', &
2317                         ' on host "', host, '"'
2318          ENDIF
2319       ENDIF
2320    ENDIF
2321
2322
2323 END SUBROUTINE check_parameters
Note: See TracBrowser for help on using the repository browser.