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

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

Fixes for spinup and sky-view-factor (SVF) treatment

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