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

Last change on this file since 2012 was 2012, checked in by kanani, 8 years ago

last commit documented

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