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

Last change on this file since 2696 was 2696, checked in by kanani, 6 years ago

Merge of branch palm4u into trunk

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