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

Last change on this file since 3313 was 3313, checked in by knoop, 5 years ago

Added test mode to palmrun

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