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

Last change on this file since 2375 was 2375, checked in by schwenkel, 7 years ago

improved aerosol initialization for bulk microphysics

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