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

Last change on this file since 2038 was 2038, checked in by knoop, 7 years ago

last commit documented

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