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

Last change on this file since 2232 was 2232, checked in by suehring, 8 years ago

Adjustments according new topography and surface-modelling concept implemented

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