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

Last change on this file since 2317 was 2310, checked in by gronemeier, 7 years ago

bugfix: re-arranged call for error messages for ENVPAR file

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