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

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

introduced new module date_and_time_mod

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