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

Last change on this file since 2575 was 2575, checked in by maronga, 6 years ago

bugfix in radiation model and improvements in land surface scheme

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