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

Last change on this file since 3083 was 3083, checked in by gronemeier, 6 years ago

merge with branch rans

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