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

Last change on this file since 3353 was 3347, checked in by suehring, 6 years ago

Offline nesting revised and separated from large_scale_forcing_mod; Time-dependent synthetic turbulence generator; bugfixes in USM/LSM radiation model and init_pegrid

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