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

Last change on this file since 2968 was 2967, checked in by raasch, 6 years ago

bugfix: missing parallel cpp-directives added

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