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

Last change on this file since 2478 was 2397, checked in by suehring, 7 years ago

Enable initialization of 3d model by user in the child domain

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