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

Last change on this file since 2881 was 2881, checked in by maronga, 6 years ago

bugfix in land surface model and new option in spinup

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