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

Last change on this file since 1955 was 1955, checked in by hellstea, 8 years ago

Bugfix on nest-domain

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