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

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

Humidity and passive scalar also separated in nesting mode

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