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

Last change on this file since 2124 was 2119, checked in by raasch, 7 years ago

last commit documented

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