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

Last change on this file since 3250 was 3248, checked in by sward, 6 years ago

Minor format changes

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