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

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