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

Last change on this file since 2716 was 2716, checked in by kanani, 6 years ago

Correction of "Former revisions" section

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