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

Last change on this file since 3987 was 3987, checked in by kanani, 5 years ago

clean up location, debug and error messages

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