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

Last change on this file since 3472 was 3472, checked in by suehring, 5 years ago

module for virtual measurements added (in a preliminary state); new public routines to input NetCDF data directly from modules

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