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

Last change on this file since 3046 was 3046, checked in by Giersch, 6 years ago

Remaining error messages revised, comments extended

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