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

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

last commit documented

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