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

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

last commit documented

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