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

Last change on this file since 2898 was 2894, checked in by Giersch, 6 years ago

Reading/Writing? data in case of restart runs revised

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