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

Last change on this file since 3198 was 3184, checked in by suehring, 6 years ago

Bugfix in setting of passive scalar and chemistry top boundary condary condition in offline nesting mode

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