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

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

Removal of chem directive, plus minor changes

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