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

Last change on this file since 3780 was 3747, checked in by gronemeier, 5 years ago

bugfix: set user-defined statistic regions within new routine user_init_arrays and initialize rmask before this routine is called. This ensures that the correct number of grid points is calculated for each statistic region within init_3d_model.
(init_3d_model, module_interface, user_module)
bugfix: update former revisions section (modules, parin)

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