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

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

Improvement for spinup checks concering consistency in nesting mode

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