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

Last change on this file since 2257 was 2233, checked in by suehring, 7 years ago

last commit documented

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