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

Last change on this file since 2000 was 2000, checked in by knoop, 8 years ago

Forced header and separation lines into 80 columns

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