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

Last change on this file since 2296 was 2296, checked in by maronga, 7 years ago

added new spinup mechanism for surface/radiation models

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