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

Last change on this file since 2918 was 2906, checked in by Giersch, 6 years ago

new procedure for reading/writing svf data, initialization of local variable ids

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