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

Last change on this file since 2799 was 2773, checked in by suehring, 6 years ago

Nesting for chemical species implemented; Bugfix passive scalar boundary condition after anterpolation; Timeseries output of surface temperature; Enable initialization of 3D topography (was commented out so far)

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