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

Last change on this file since 2365 was 2365, checked in by kanani, 7 years ago

Vertical nesting implemented (SadiqHuq?)

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