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

Last change on this file since 2550 was 2550, checked in by boeske, 7 years ago

enable simulations with complex terrain

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