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

Last change on this file since 2826 was 2826, checked in by hellstea, 6 years ago

Bugfix in setting the default boundary conditions for nest domains

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