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

Last change on this file since 1956 was 1956, checked in by hellstea, 8 years ago

last commit documented

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