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

Last change on this file since 3232 was 3204, checked in by raasch, 6 years ago

additional check for nz

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