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

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

corrected timestamp in header

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