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

Last change on this file since 2298 was 2298, checked in by raasch, 7 years ago

write_binary is of type LOGICAL now, MPI2-related code removed, obsolete variables removed, sendrecv_in_background related parts removed, missing variable descriptions added

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