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

Last change on this file since 2372 was 2372, checked in by sward, 7 years ago

y_shift for periodic boundary conditions

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