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

Last change on this file since 3274 was 3274, checked in by knoop, 6 years ago

Modularization of all bulk cloud physics code components

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