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

Last change on this file since 3646 was 3646, checked in by kanani, 5 years ago

Bugfix: replace simulated_time by time_since_reference_point where required

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