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

Last change on this file since 2563 was 2563, checked in by Giersch, 7 years ago

Restart runs with the usage of the wind turbine model are possible now. Further small at reading/writing restart data

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