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

Last change on this file since 2037 was 2037, checked in by knoop, 7 years ago

Anelastic approximation implemented

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