source: palm/trunk/SOURCE/parin.f90 @ 2263

Last change on this file since 2263 was 2259, checked in by gronemeier, 7 years ago

Implemented synthetic turbulence generator

  • Property svn:keywords set to Id
File size: 27.4 KB
Line 
1!> @file parin.f90
2!------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2017 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: parin.f90 2259 2017-06-08 09:09:11Z schwenkel $
27! Implemented synthetic turbulence generator
28!
29! 2233 2017-05-30 18:08:54Z suehring
30!
31! 2232 2017-05-30 17:47:52Z suehring
32! typo corrected
33! +wall_salinityflux
34! +tunnel_height, tunnel_lenght, tunnel_width_x, tunnel_width_y,
35!  tunnel_wall_depth
36!
37! 2118 2017-01-17 16:38:49Z raasch
38! -background_communication from inipar
39!
40! 2050 2016-11-08 15:00:55Z gronemeier
41! Implement turbulent outflow condition
42!
43! 2037 2016-10-26 11:15:40Z knoop
44! Anelastic approximation implemented
45!
46! 2035 2016-10-24 15:06:17Z suehring
47! Remove check for npex and npey in nesting case
48!
49! 2011 2016-09-19 17:29:57Z kanani
50! Added flag lsf_exception to allow explicit enabling of large scale forcing
51! together with buildings on flat terrain.
52!
53! 2007 2016-08-24 15:47:17Z kanani
54! Added call to urban surface model for reading of &urban_surface_par
55!
56! 2004 2016-08-24 10:25:59Z suehring
57! Humidity and passive scalar treated separately in nesting mode
58!
59! 2000 2016-08-20 18:09:15Z knoop
60! Forced header and separation lines into 80 columns
61!
62! 1992 2016-08-12 15:14:59Z suehring
63! +top_scalarflux
64!
65! 1960 2016-07-12 16:34:24Z suehring
66! Allocation of s_init
67!
68! 1957 2016-07-07 10:43:48Z suehring
69! flight module added
70!
71! 1955 2016-07-01 12:38:59Z hellstea
72! The parameter intializating_actions is set to 'set_constant_profiles for
73! all nest domains in order to make sure that diagnostic variables are properly
74! initialized for nest domains. Prognostic variables are later initialized by
75! interpolation from the parent domain.
76!
77! 1917 2016-05-27 14:28:12Z witha
78! Initial version of purely vertical nesting introduced.
79!
80! 1914 2016-05-26 14:44:07Z witha
81! Added call to wind turbine model for reading of &wind_turbine_par
82!
83! 1849 2016-04-08 11:33:18Z hoffmann
84! Adapted for modularization of microphysics
85!
86! 1833 2016-04-07 14:23:03Z raasch
87! call of spectra_parin
88!
89! 1831 2016-04-07 13:15:51Z hoffmann
90! turbulence renamed collision_turbulence, drizzle renamed
91! cloud_water_sedimentation
92! curvature_solution_effects removed
93!
94! 1826 2016-04-07 12:01:39Z maronga
95! Added call to radiation model for reading of &radiation_par.
96! Added call to plant canopy model for reading of &canopy_par.
97!
98! 1817 2016-04-06 15:44:20Z maronga
99! Added call to land surface model for reading of &lsm_par
100!
101! 1804 2016-04-05 16:30:18Z maronga
102! Removed code for parameter file check (__check)
103!
104! 1783 2016-03-06 18:36:17Z raasch
105! +netcdf_deflate in d3par, netcdf module and variable names changed
106!
107! 1764 2016-02-28 12:45:19Z raasch
108! cpp-statements for nesting removed, explicit settings of boundary conditions
109! in nest domains,
110! bugfix: npex/npey message moved from inipar to d3par
111! bugfix: check of lateral boundary conditions from check_parameters to here,
112! because they will be already used in init_pegrid and init_grid
113!
114! 1762 2016-02-25 12:31:13Z hellstea
115! Introduction of nested domain feature
116!
117! 1691 2015-10-26 16:17:44Z maronga
118! Added parameter most_method. Renamed prandtl_layer to constant_flux_layer.
119!
120! 1682 2015-10-07 23:56:08Z knoop
121! Code annotations made doxygen readable
122!
123! 1560 2015-03-06 10:48:54Z keck
124! +recycling_yshift
125!
126! 1496 2014-12-02 17:25:50Z maronga
127! Renamed: "radiation -> "cloud_top_radiation"
128!
129! 1484 2014-10-21 10:53:05Z kanani
130! Changes due to new module structure of the plant canopy model:
131!   canopy-model related parameters moved to new package canopy_par in
132!   subroutine package_parin
133!
134! 1429 2014-07-15 12:53:45Z knoop
135! +ensemble_member_nr to prepare the random_generator for ensemble runs
136!
137! 1402 2014-05-09 14:25:13Z raasch
138! location messages modified, batch_job included in envpar-NAMELIST
139!
140! 1384 2014-05-02 14:31:06Z raasch
141! location messages added
142!
143! 1365 2014-04-22 15:03:56Z boeske
144! Usage of large scale forcing enabled:
145! +use_subsidence_tendencies
146!
147! 1361 2014-04-16 15:17:48Z hoffmann
148! +call_microphysics_at_all_substeps
149!
150! 1359 2014-04-11 17:15:14Z hoffmann
151! REAL constants provided with KIND-attribute
152!
153! 1353 2014-04-08 15:21:23Z heinze
154! REAL constants provided with KIND-attribute
155!
156! 1327 2014-03-21 11:00:16Z raasch
157! -data_output_format, do3d_compress, do3d_comp_prec
158!
159! 1320 2014-03-20 08:40:49Z raasch
160! ONLY-attribute added to USE-statements,
161! kind-parameters added to all INTEGER and REAL declaration statements,
162! kinds are defined in new module kinds,
163! old module precision_kind is removed,
164! revision history before 2012 removed,
165! comment fields (!:) to be used for variable explanations added to
166! all variable declaration statements
167!
168! 1318 2014-03-17 13:35:16Z raasch
169! +cpu_log_barrierwait in d3par
170!
171! 1301 2014-03-06 13:29:46Z heinze
172! +large_scale_subsidence
173!
174! 1241 2013-10-30 11:36:58Z heinze
175! +nudging
176! +large_scale_forcing
177!
178! 1216 2013-08-26 09:31:42Z raasch
179! +transpose_compute_overlap in inipar
180!
181! 1195 2013-07-01 12:27:57Z heinze
182! Bugfix: allocate ref_state
183!
184! 1179 2013-06-14 05:57:58Z raasch
185! +reference_state in inipar
186!
187! 1159 2013-05-21 11:58:22Z fricke
188! +use_cmax
189!
190! 1128 2013-04-12 06:19:32Z raasch
191! +background_communication in inipar
192!
193! 1115 2013-03-26 18:16:16Z hoffmann
194! unused variables removed
195!
196! 1092 2013-02-02 11:24:22Z raasch
197! unused variables removed
198!
199! 1065 2012-11-22 17:42:36Z hoffmann
200! +nc, c_sedimentation, limiter_sedimentation, turbulence
201! -mu_constant, mu_constant_value
202!
203! 1053 2012-11-13 17:11:03Z hoffmann
204! necessary expansions according to the two new prognostic equations (nr, qr)
205! of the two-moment cloud physics scheme and steering parameters:
206! +*_init, *_surface, *_surface_initial_change, *_vertical_gradient,
207! +*_vertical_gradient_level, surface_waterflux_*,
208! +cloud_scheme, drizzle, mu_constant, mu_constant_value, ventilation_effect
209!
210! 1036 2012-10-22 13:43:42Z raasch
211! code put under GPL (PALM 3.9)
212!
213! 1015 2012-09-27 09:23:24Z raasch
214! -adjust_mixing_length
215!
216! 1003 2012-09-14 14:35:53Z raasch
217! -grid_matching
218!
219! 1001 2012-09-13 14:08:46Z raasch
220! -cut_spline_overshoot, long_filter_factor, overshoot_limit_*, ups_limit_*
221!
222! 996 2012-09-07 10:41:47Z raasch
223! -use_prior_plot1d_parameters
224!
225! 978 2012-08-09 08:28:32Z fricke
226! -km_damp_max, outflow_damping_width
227! +pt_damping_factor, pt_damping_width
228! +z0h_factor
229!
230! 964 2012-07-26 09:14:24Z raasch
231! -cross_normalized_x, cross_normalized_y, cross_xtext, z_max_do1d,
232! z_max_do1d_normalized
233!
234! 940 2012-07-09 14:31:00Z raasch
235! +neutral in inipar
236!
237! 927 2012-06-06 19:15:04Z raasch
238! +masking_method in inipar
239!
240! 824 2012-02-17 09:09:57Z raasch
241! +curvature_solution_effects in inipar
242!
243! 809 2012-01-30 13:32:58Z maronga
244! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
245!
246! 807 2012-01-25 11:53:51Z maronga
247! New cpp directive "__check" implemented which is used by check_namelist_files
248!
249! Revision 1.1  1997/07/24 11:22:50  raasch
250! Initial revision
251!
252!
253! Description:
254! ------------
255!> This subroutine reads variables controling the run from the NAMELIST files
256!------------------------------------------------------------------------------!
257 SUBROUTINE parin
258 
259
260    USE arrays_3d,                                                             &
261        ONLY:  pt_init, q_init, ref_state, s_init, sa_init, ug, u_init, v_init,&
262               vg
263
264    USE plant_canopy_model_mod,                                                &
265         ONLY: pcm_parin 
266
267    USE control_parameters
268
269    USE cpulog,                                                                &
270        ONLY:  cpu_log_barrierwait
271
272    USE dvrp_variables,                                                        &
273        ONLY:  local_dvrserver_running
274
275    USE flight_mod,                                                            &
276        ONLY:  flight_parin
277
278    USE grid_variables,                                                        &
279        ONLY:  dx, dy
280
281    USE indices,                                                               &
282        ONLY:  nx, ny, nz
283
284    USE kinds
285
286    USE land_surface_model_mod,                                                &
287        ONLY: lsm_parin
288
289    USE microphysics_mod,                                                      &
290        ONLY:  c_sedimentation, cloud_water_sedimentation,                     &
291               collision_turbulence, limiter_sedimentation, nc_const,          &
292               ventilation_effect
293
294    USE model_1d,                                                              &
295        ONLY:  damp_level_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
296
297    USE pegrid
298
299    USE netcdf_interface,                                                      &
300        ONLY:  netcdf_data_format, netcdf_deflate, netcdf_precision
301
302    USE pmc_interface,                                                         &
303        ONLY:  nested_run, nesting_mode
304
305    USE profil_parameter,                                                      &
306        ONLY:  cross_profiles, cross_ts_uymax, cross_ts_uymin,                 &
307               profile_columns, profile_rows
308
309    USE progress_bar,                                                          &
310        ONLY :  batch_job
311
312    USE radiation_model_mod,                                                   &
313        ONLY: radiation_parin 
314
315    USE spectra_mod,                                                           &
316        ONLY :  spectra_parin
317
318    USE statistics,                                                            &
319        ONLY:  hom, hom_sum, pr_palm, region, statistic_regions
320
321    USE synthetic_turbulence_generator_mod,                                    &
322        ONLY:  stg_parin
323
324    USE urban_surface_mod,                                                     &
325        ONLY: usm_parin
326
327    USE wind_turbine_model_mod,                                                &
328        ONLY:  wtm_parin
329
330
331    IMPLICIT NONE
332
333    INTEGER(iwp) ::  i   !<
334
335
336    NAMELIST /inipar/  alpha_surface, approximation, bc_e_b, bc_lr,            &
337                       bc_ns, bc_p_b, bc_p_t, bc_pt_b, bc_pt_t, bc_q_b,        &
338             bc_q_t,bc_s_b, bc_s_t, bc_sa_t, bc_uv_b, bc_uv_t,                 &
339             bottom_salinityflux, building_height, building_length_x,          &
340             building_length_y, building_wall_left, building_wall_south,       &
341             call_psolver_at_all_substeps, call_microphysics_at_all_substeps,  &
342             canyon_height,                                                    &
343             canyon_width_x, canyon_width_y, canyon_wall_left,                 &
344             canyon_wall_south, c_sedimentation, cfl_factor, cloud_droplets,   &
345             cloud_physics, cloud_scheme, cloud_top_radiation,                 &
346             cloud_water_sedimentation,                                        &
347             collective_wait, collision_turbulence, conserve_volume_flow,      &
348             conserve_volume_flow_mode, constant_flux_layer,                   &
349             coupling_start_time,                                              &
350             cycle_mg, damp_level_1d,                                          &
351             dissipation_1d,                                                   &
352             dp_external, dp_level_b, dp_smooth, dpdxy,                        &
353             dt, dt_pr_1d, dt_run_control_1d, dx, dy, dz, dz_max,              &
354             dz_stretch_factor, dz_stretch_level, end_time_1d,                 &
355             ensemble_member_nr, e_init, e_min, fft_method,                    &
356             flux_input_mode, flux_output_mode,                                &
357             galilei_transformation, humidity,                                 &
358             inflow_damping_height, inflow_damping_width,                      &
359             inflow_disturbance_begin, inflow_disturbance_end,                 &
360             initializing_actions, km_constant,                                &
361             large_scale_forcing, large_scale_subsidence,                      &
362             limiter_sedimentation,                                            &
363             loop_optimization, lsf_exception, masking_method, mg_cycles,      &
364             mg_switch_to_pe0_level, mixing_length_1d, momentum_advec,         &
365             most_method, nc_const, netcdf_precision, neutral, ngsrb,          &
366             nsor, nsor_ini, nudging, nx, ny, nz, ocean, omega, omega_sor,     &
367             outflow_source_plane, passive_scalar, phi,                        &
368             prandtl_number, precipitation, psolver, pt_damping_factor,        &
369             pt_damping_width, pt_reference, pt_surface,                       &
370             pt_surface_initial_change, pt_vertical_gradient,                  &
371             pt_vertical_gradient_level, q_surface, q_surface_initial_change,  &
372             q_vertical_gradient, q_vertical_gradient_level,                   &
373             random_generator, random_heatflux,                                &
374             rayleigh_damping_factor, rayleigh_damping_height,                 &
375             recycling_width, recycling_yshift,                                &
376             reference_state, residual_limit,                                  &
377             roughness_length, sa_surface,                                     &
378             sa_vertical_gradient, sa_vertical_gradient_level, scalar_advec,   &
379             scalar_rayleigh_damping,                                          &
380             statistic_regions, subs_vertical_gradient,                        &
381             subs_vertical_gradient_level, surface_heatflux, surface_pressure, &
382             surface_scalarflux, surface_waterflux,                            &
383             s_surface, s_surface_initial_change, s_vertical_gradient,         &
384             s_vertical_gradient_level, timestep_scheme,                       &
385             topography, topography_grid_convention, top_heatflux,             &
386             top_momentumflux_u, top_momentumflux_v, top_salinityflux,         &
387             top_scalarflux, transpose_compute_overlap,                        &
388             tunnel_height, tunnel_length, tunnel_width_x, tunnel_width_y,     &
389             tunnel_wall_depth,                                                &
390             turbulent_inflow, turbulent_outflow,                              &
391             use_subsidence_tendencies, ug_surface, ug_vertical_gradient,      &
392             ug_vertical_gradient_level, use_surface_fluxes, use_cmax,         &
393             use_top_fluxes, use_ug_for_galilei_tr, use_upstream_for_tke,      &
394             uv_heights, u_bulk, u_profile, vg_surface, vg_vertical_gradient,  &
395             vg_vertical_gradient_level, v_bulk, v_profile, ventilation_effect,&
396             wall_adjustment, wall_heatflux, wall_humidityflux,                &
397             wall_salinityflux, wall_scalarflux, zeta_max, zeta_min, z0h_factor
398     
399    NAMELIST /d3par/  averaging_interval, averaging_interval_pr,               &
400             cpu_log_barrierwait, create_disturbances,                         &
401             cross_profiles, cross_ts_uymax, cross_ts_uymin,                   &
402             data_output, data_output_masks,                                   &
403             data_output_pr, data_output_2d_on_each_pe, disturbance_amplitude, &
404             disturbance_energy_limit, disturbance_level_b,                    &
405             disturbance_level_t, do2d_at_begin, do3d_at_begin,                &
406             dt, dt_averaging_input, dt_averaging_input_pr,                    &
407             dt_coupling, dt_data_output, dt_data_output_av, dt_disturb,       &
408             dt_domask, dt_dopr, dt_dopr_listing, dt_dots, dt_do2d_xy,         &
409             dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_max, dt_restart,              &
410             dt_run_control,end_time, force_print_header, mask_scale_x,        &
411             mask_scale_y, mask_scale_z, mask_x, mask_y, mask_z, mask_x_loop,  &
412             mask_y_loop, mask_z_loop, netcdf_data_format, netcdf_deflate,     &
413             normalizing_region, npex, npey, nz_do3d,                          &
414             precipitation_amount_interval, profile_columns, profile_rows,     &
415             restart_time, section_xy, section_xz, section_yz,                 &
416             skip_time_data_output, skip_time_data_output_av, skip_time_dopr,  &
417             skip_time_do2d_xy, skip_time_do2d_xz, skip_time_do2d_yz,          &
418             skip_time_do3d, skip_time_domask, synchronous_exchange,           &
419             termination_time_needed, z_max_do2d
420
421
422    NAMELIST /envpar/  batch_job, host, local_dvrserver_running,               &
423                       maximum_cpu_time_allowed, maximum_parallel_io_streams,  &
424                       revision, return_addres, return_username,               &
425                       run_identifier, tasks_per_node, write_binary
426
427!
428!-- First read values of environment variables (this NAMELIST file is
429!-- generated by mrun)
430    CALL location_message( 'reading environment parameters from ENVPAR', .FALSE. )
431    OPEN ( 90, FILE='ENVPAR', STATUS='OLD', FORM='FORMATTED', ERR=30 )
432    READ ( 90, envpar, ERR=31, END=32 )
433    CLOSE ( 90 )
434    CALL location_message( 'finished', .TRUE. )
435
436!
437!-- Calculate the number of groups into which parallel I/O is split.
438!-- The default for files which are opened by all PEs (or where each
439!-- PE opens his own independent file) is, that all PEs are doing input/output
440!-- in parallel at the same time. This might cause performance or even more
441!-- severe problems depending on the configuration of the underlying file
442!-- system.
443!-- First, set the default:
444    IF ( maximum_parallel_io_streams == -1  .OR.                               &
445         maximum_parallel_io_streams > numprocs )  THEN
446       maximum_parallel_io_streams = numprocs
447    ENDIF
448!
449!-- Now calculate the number of io_blocks and the io_group to which the
450!-- respective PE belongs. I/O of the groups is done in serial, but in parallel
451!-- for all PEs belonging to the same group.
452!-- These settings are repeated in init_pegrid for the communicator comm2d,
453!-- which is not available here
454    io_blocks = numprocs / maximum_parallel_io_streams
455    io_group  = MOD( myid+1, io_blocks )
456
457    CALL location_message( 'reading NAMELIST parameters from PARIN', .FALSE. )
458!
459!-- Data is read in parallel by groups of PEs
460    DO  i = 0, io_blocks-1
461       IF ( i == io_group )  THEN
462
463!
464!--       Open the NAMELIST-file which is send with this job
465          CALL check_open( 11 )
466
467!
468!--       Read the control parameters for initialization.
469!--       The namelist "inipar" must be provided in the NAMELIST-file.
470          READ ( 11, inipar, ERR=10, END=11 )
471
472          GOTO 12
473
474 10       message_string = 'errors in \$inipar &or no \$inipar-namelist ' //   &
475                           'found (CRAY-machines only)'
476          CALL message( 'parin', 'PA0271', 1, 2, 0, 6, 0 )
477
478 11       message_string = 'no \$inipar-namelist found'
479          CALL message( 'parin', 'PA0272', 1, 2, 0, 6, 0 )
480
481!
482!--       If required, read control parameters from restart file (produced by
483!--       a prior run). All PEs are reading from file created by PE0 (see
484!--       check_open)
485 12       IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
486             CALL read_var_list
487!
488!--          The restart file will be reopened when reading the subdomain data
489             CALL close_file( 13 )
490
491!
492!--          Increment the run count
493             runnr = runnr + 1
494          ENDIF
495
496!
497!--       In case of nested runs, explicitly set nesting boundary conditions.
498!--       This will overwrite the user settings. bc_lr and bc_ns always need
499!--       to be cyclic for vertical nesting.
500          IF ( nesting_mode == 'vertical' )  THEN
501             IF (bc_lr /= 'cyclic' .OR. bc_ns /= 'cyclic' )  THEN
502                WRITE ( message_string, *)  'bc_lr and bc_ns were set to ,',   &
503                                      'cyclic for vertical nesting'
504                CALL message( 'parin', 'PA0428', 0, 0, 0, 6, 0 )
505                bc_lr   = 'cyclic'
506                bc_ns   = 'cyclic'
507             ENDIF
508             IF ( nest_domain )  THEN
509                bc_uv_t = 'nested'
510                bc_pt_t = 'nested'
511                bc_q_t  = 'nested'
512                bc_s_t  = 'nested'
513                bc_p_t  = 'neumann'
514             ENDIF
515          ELSE
516         
517!
518!--       For other nesting modes only set boundary conditions for
519!--       nested domains.
520             IF ( nest_domain )  THEN
521                bc_lr   = 'nested'
522                bc_ns   = 'nested'
523                bc_uv_t = 'nested'
524                bc_pt_t = 'nested'
525                bc_q_t  = 'nested'
526                bc_s_t  = 'nested'
527                bc_p_t  = 'neumann'
528             ENDIF
529          ENDIF
530
531!         
532!--       In case of nested runs, make sure that initializing_actions =
533!--       'set_constant_profiles' even though the constant-profiles
534!--       initializations for the prognostic variables will be overwritten
535!--       by pmci_child_initialize and pmci_parent_initialize. This is,
536!--       however, important e.g. to make sure that diagnostic variables
537!--       are set properly.
538          IF ( nest_domain )  THEN
539             initializing_actions = 'set_constant_profiles'
540          ENDIF
541
542           
543!
544!--       Check validity of lateral boundary conditions. This has to be done
545!--       here because they are already used in init_pegrid and init_grid and
546!--       therefore cannot be check in check_parameters
547          IF ( bc_lr /= 'cyclic'  .AND.  bc_lr /= 'dirichlet/radiation'  .AND. &
548               bc_lr /= 'radiation/dirichlet'  .AND.  bc_lr /= 'nested' )  THEN
549             message_string = 'unknown boundary condition: bc_lr = "' // &
550                              TRIM( bc_lr ) // '"'
551             CALL message( 'check_parameters', 'PA0049', 1, 2, 0, 6, 0 )
552          ENDIF
553          IF ( bc_ns /= 'cyclic'  .AND.  bc_ns /= 'dirichlet/radiation'  .AND. &
554               bc_ns /= 'radiation/dirichlet'  .AND.  bc_ns /= 'nested' )  THEN
555             message_string = 'unknown boundary condition: bc_ns = "' // &
556                              TRIM( bc_ns ) // '"'
557             CALL message( 'check_parameters', 'PA0050', 1, 2, 0, 6, 0 )
558          ENDIF
559
560!
561!--       Set internal variables used for speed optimization in if clauses
562          IF ( bc_lr /= 'cyclic' )               bc_lr_cyc    = .FALSE.
563          IF ( bc_lr == 'dirichlet/radiation' )  bc_lr_dirrad = .TRUE.
564          IF ( bc_lr == 'radiation/dirichlet' )  bc_lr_raddir = .TRUE.
565          IF ( bc_ns /= 'cyclic' )               bc_ns_cyc    = .FALSE.
566          IF ( bc_ns == 'dirichlet/radiation' )  bc_ns_dirrad = .TRUE.
567          IF ( bc_ns == 'radiation/dirichlet' )  bc_ns_raddir = .TRUE.
568
569!
570!--       Definition of names of areas used for computing statistics. They must
571!--       be defined at this place, because they are allowed to be redefined by
572!--       the user in user_parin.
573          region = 'total domain'
574
575!
576!--       Read runtime parameters given by the user for this run (namelist
577!--       "d3par"). The namelist "d3par" can be omitted. In that case, default
578!--       values are used for the parameters.
579          READ ( 11, d3par, END=20 )
580
581!
582!--       Check if land surface model is used and read &lsm_par if required
583          CALL lsm_parin
584
585!
586!--       Check if urban surface model is used and read &urban_surface_par if required
587          CALL usm_parin
588
589!
590!--       Check if spectra shall be calculated and read spectra_par if required
591          CALL spectra_parin
592
593!
594!--       Check if radiation model is used and read &radiation_par if required
595          CALL radiation_parin
596 
597 
598!--       Check if plant canopy model is used and read &canopy_par if required
599          CALL pcm_parin
600 
601!
602!--       Read control parameters for optionally used model software packages
603 20       CALL package_parin
604
605!
606!--       Check if wind turbine model is used and read &wind_turbine_par if
607!--       required
608          CALL wtm_parin
609!
610!--       Check if virtual flights should be carried out and read &flight_par
611!--       if required
612          CALL flight_parin
613!
614!--       Check if synthetic turbulence generator is used and read stg_par if
615!--       required
616          CALL stg_parin
617!
618!--       Read user-defined variables
619          CALL user_parin
620
621!
622!--       Check in case of initial run, if the grid point numbers are well
623!--       defined and allocate some arrays which are already needed in
624!--       init_pegrid or check_parameters. During restart jobs, these arrays
625!--       will be allocated in read_var_list. All other arrays are allocated
626!--       in init_3d_model.
627          IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
628
629             IF ( nx <= 0 )  THEN
630                WRITE( message_string, * ) 'no value or wrong value given',    &
631                                           ' for nx: nx=', nx
632                CALL message( 'parin', 'PA0273', 1, 2, 0, 6, 0 )
633             ENDIF
634             IF ( ny <= 0 )  THEN
635                WRITE( message_string, * ) 'no value or wrong value given',    &
636                                           ' for ny: ny=', ny
637                CALL message( 'parin', 'PA0274', 1, 2, 0, 6, 0 )
638             ENDIF
639             IF ( nz <= 0 )  THEN
640                WRITE( message_string, * ) 'no value or wrong value given',    &
641                                           ' for nz: nz=', nz
642                CALL message( 'parin', 'PA0275', 1, 2, 0, 6, 0 )
643             ENDIF
644!
645!--          ATTENTION: in case of changes to the following statement please
646!--                  also check the allocate statement in routine read_var_list
647             ALLOCATE( pt_init(0:nz+1), q_init(0:nz+1), s_init(0:nz+1),        &
648                       ref_state(0:nz+1), sa_init(0:nz+1), ug(0:nz+1),         &
649                       u_init(0:nz+1), v_init(0:nz+1), vg(0:nz+1),             &
650                       hom(0:nz+1,2,pr_palm+max_pr_user,0:statistic_regions),  &
651                       hom_sum(0:nz+1,pr_palm+max_pr_user,0:statistic_regions) )
652
653             hom = 0.0_wp
654
655          ENDIF
656
657!
658!--       NAMELIST-file is not needed anymore
659          CALL close_file( 11 )
660
661       ENDIF
662#if defined( __parallel )
663       CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
664#endif
665    ENDDO
666
667    CALL location_message( 'finished', .TRUE. )
668
669    RETURN
670
671 30 message_string = 'local file ENVPAR not found' //                          &
672                     '&some variables for steering may not be properly set'
673    CALL message( 'parin', 'PA0276', 0, 1, 0, 6, 0 )
674    RETURN
675
676 31 message_string = 'errors in local file ENVPAR' //                          &
677                     '&some variables for steering may not be properly set'
678    CALL message( 'parin', 'PA0277', 0, 1, 0, 6, 0 )
679    RETURN
680
681 32 message_string = 'no envpar-NAMELIST found in local file ENVPAR'  //       &
682                     '&some variables for steering may not be properly set'
683    CALL message( 'parin', 'PA0278', 0, 1, 0, 6, 0 )
684
685 END SUBROUTINE parin
Note: See TracBrowser for help on using the repository browser.