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

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

Bugfix: remove skipping of namelists in case of omitted d3par

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