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

Last change on this file since 3562 was 3545, checked in by gronemeier, 5 years ago

removed debug output; removed rans_mode from namelist; altered order of check_parameter-calls of modules

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