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

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

last commit documented

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