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

Last change on this file since 2067 was 2051, checked in by gronemeier, 7 years ago

last commit documented

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