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

Last change on this file since 2938 was 2938, checked in by suehring, 6 years ago

Nesting in RANS-LES and RANS-RANS mode enabled; synthetic turbulence generator at all lateral boundaries in nesting or non-cyclic forcing mode; revised Inifor initialization in nesting mode

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