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

Last change on this file since 2008 was 2008, checked in by kanani, 8 years ago

last commit documented

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