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

Last change on this file since 2823 was 2817, checked in by knoop, 6 years ago

Preliminary gust module interface implemented

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