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

Last change on this file since 3435 was 3435, checked in by gronemeier, 6 years ago

new: terrain-following masked output; bugfixes: increase vertical dimension of gamma_w_green_sat by 1, add checks for masked output for chemistry_model and radiation_model, reordered calls to xxx_define_netcdf_grid in masked output part

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