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

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

renamed all Fortran NAMELISTS

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