source: palm/trunk/SOURCE/init_3d_model.f90 @ 3288

Last change on this file since 3288 was 3288, checked in by suehring, 6 years ago

Bugfix in allocation of mean_inflow_profiles in case of restarts

  • Property svn:keywords set to Id
  • Property svn:mergeinfo set to False
    /palm/branches/forwind/SOURCE/init_3d_model.f901564-1913
    /palm/branches/palm4u/SOURCE/init_3d_model.f902540-2692
    /palm/branches/rans/SOURCE/init_3d_model.f902078-3128
File size: 96.1 KB
RevLine 
[1682]1!> @file init_3d_model.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!
[254]20! Current revisions:
[732]21! ------------------
[3288]22! Introduce module parameter for number of inflow profiles
[2233]23!
24! Former revisions:
25! -----------------
26! $Id: init_3d_model.f90 3288 2018-09-28 10:23:08Z suehring $
[3274]27! Modularization of all bulk cloud physics code components
28!
29! 3241 2018-09-12 15:02:00Z raasch
[3241]30! unused variables removed
31!
32! 3234 2018-09-07 13:46:58Z schwenkel
[3234]33! The increase of dots_num in case of radiation or land surface model must
34! be done before user_init is called
35!
36! 3183 2018-07-27 14:25:55Z suehring
[3183]37! Revise Inifor initialization
38!
39! 3182 2018-07-27 13:36:03Z suehring
[3159]40! Added multi agent system
41!
42! 3129 2018-07-16 07:45:13Z gronemeier
[3051]43! Move initialization call for nudging and 1D/3D offline nesting.
44! Revise initialization with inifor data.
45!
46! 3045 2018-05-28 07:55:41Z Giersch
[3045]47! Error messages revised
48!
[3049]49! 3045 2018-05-28 07:55:41Z Giersch
50! Error messages revised
51!
[3045]52! 3042 2018-05-25 10:44:37Z schwenkel
[3042]53! Changed the name specific humidity to mixing ratio
54!
55! 3040 2018-05-25 10:22:08Z schwenkel
[3035]56! Add option to initialize warm air bubble close to surface
57!
58! 3014 2018-05-09 08:42:38Z maronga
[3014]59! Bugfix: initialization of ts_value missing
60!
61! 3011 2018-05-07 14:38:42Z schwenkel
[3011]62! removed redundant if statement
63!
64! 3004 2018-04-27 12:33:25Z Giersch
[3004]65! precipitation_rate removed
66!
67! 2995 2018-04-19 12:13:16Z Giersch
[2995]68! CALL radiation_control is not necessary during initialization because
69! calculation of radiative fluxes at model start is done in radiation_init
70! in any case
71!
72! 2977 2018-04-17 10:27:57Z kanani
[2977]73! Implement changes from branch radiation (r2948-2971) with minor modifications
74! (moh.hefny):
75! - set radiation_interactions according to the existence of urban/land vertical
76!   surfaces and trees to activiate RTM
77! - set average_radiation to TRUE if RTM is activiated
78!
79! 2938 2018-03-27 15:52:42Z suehring
[2938]80! - Revise Inifor initialization for geostrophic wind components
81! - Initialize synthetic turbulence generator in case of Inifor initialization 
82!
83! 2936 2018-03-27 14:49:27Z suehring
[2934]84! Synchronize parent and child models after initialization.
85! Remove obsolete masking of topography grid points for Runge-Kutta weighted
86! tendency arrays.
87!
88! 2920 2018-03-22 11:22:01Z kanani
[2920]89! Add call for precalculating apparent solar positions (moh.hefny)
90!
91! 2906 2018-03-19 08:56:40Z Giersch
[2906]92! The variables read/write_svf_on_init have been removed. Instead ENVIRONMENT
93! variables read/write_svf have been introduced. Location_message has been
94! added.
95!
96! 2894 2018-03-15 09:17:58Z Giersch
[2894]97! Renamed routines with respect to reading restart data, file 13 is closed in
98! rrd_read_parts_of_global now
99!
100! 2867 2018-03-09 09:40:23Z suehring
[2867]101! Further bugfix concerning call of user_init.
102!
103! 2864 2018-03-08 11:57:45Z suehring
[2864]104! Bugfix, move call of user_init in front of initialization of grid-point
105! arrays
106!
107! 2817 2018-02-19 16:32:21Z knoop
[2817]108! Preliminary gust module interface implemented
109!
110! 2776 2018-01-31 10:44:42Z Giersch
[2776]111! Variable use_synthetic_turbulence_generator has been abbreviated
112!
113! 2766 2018-01-22 17:17:47Z kanani
[2766]114! Removed preprocessor directive __chem
115!
116! 2758 2018-01-17 12:55:21Z suehring
[2758]117! In case of spinup of land- and urban-surface model, do not mask wind velocity
118! at first computational grid level
119!
120! 2746 2018-01-15 12:06:04Z suehring
[2746]121! Move flag plant canopy to modules
122!
123! 2718 2018-01-02 08:49:38Z maronga
[2716]124! Corrected "Former revisions" section
125!
126! 2705 2017-12-18 11:26:23Z maronga
[2705]127! Bugfix for reading initial profiles from ls/nuding file
[2716]128!
129! 2701 2017-12-15 15:40:50Z suehring
130! Changes from last commit documented
[2705]131!
[2716]132! 2700 2017-12-15 14:12:35Z suehring
[2701]133! Bugfix, missing initialization of surface attributes in case of
134! inifor-initialization branch
[2716]135!
136! 2698 2017-12-14 18:46:24Z suehring
137! Bugfix in get_topography_top_index
138!
139! 2696 2017-12-14 17:12:51Z kanani
140! Change in file header (GPL part)
[2696]141! Implementation of uv exposure model (FK)
142! Moved initialisation of diss, e, kh, km to turbulence_closure_mod (TG)
143! Added chemical emissions (FK)
144! Initialize masking arrays and number-of-grid-points arrays before initialize
145! LSM, USM and radiation module
146! Initialization with inifor (MS)
147!
148! 2618 2017-11-16 15:37:30Z suehring
[2618]149! Reorder calls of init_surfaces.
150!
151! 2564 2017-10-19 15:56:56Z Giersch
[2564]152! Variable wind_turbine was added to control_parameters.
153!
154! 2550 2017-10-16 17:12:01Z boeske
[2550]155! Modifications to cyclic fill method and turbulence recycling method in case of
156! complex terrain simulations
157!
158! 2513 2017-10-04 09:24:39Z kanani
[2513]159! Bugfix in storing initial scalar profile (wrong index)
160!
161! 2350 2017-08-15 11:48:26Z kanani
[2350]162! Bugfix in nopointer version
163!
164! 2339 2017-08-07 13:55:26Z gronemeier
[2339]165! corrected timestamp in header
166!
167! 2338 2017-08-07 12:15:38Z gronemeier
[2338]168! Modularize 1D model
169!
[2339]170! 2329 2017-08-03 14:24:56Z knoop
[2329]171! Removed temporary bugfix (r2327) as bug is properly resolved by this revision
172!
173! 2327 2017-08-02 07:40:57Z maronga
[2327]174! Temporary bugfix
175!
176! 2320 2017-07-21 12:47:43Z suehring
[2320]177! Modularize large-scale forcing and nudging
178!
179! 2292 2017-06-20 09:51:42Z schwenkel
[2292]180! Implementation of new microphysic scheme: cloud_scheme = 'morrison'
181! includes two more prognostic equations for cloud drop concentration (nc) 
182! and cloud water content (qc).
183!
184! 2277 2017-06-12 10:47:51Z kanani
[2277]185! Removed unused variable sums_up_fraction_l
186!
187! 2270 2017-06-09 12:18:47Z maronga
[2270]188! dots_num must be increased when LSM and/or radiation is used
189!
190! 2259 2017-06-08 09:09:11Z gronemeier
[2259]191! Implemented synthetic turbulence generator
192!
193! 2252 2017-06-07 09:35:37Z knoop
[2252]194! rho_air now depending on surface_pressure even in Boussinesq mode
195!
196! 2233 2017-05-30 18:08:54Z suehring
[2233]197!
198! 2232 2017-05-30 17:47:52Z suehring
[2232]199! Adjustments to new topography and surface concept:
200!   - Modify passed parameters for disturb_field
201!   - Topography representation via flags
202!   - Remove unused arrays.
203!   - Move initialization of surface-related quantities to surface_mod
[1961]204!
[2173]205! 2172 2017-03-08 15:55:25Z knoop
206! Bugfix: moved parallel random generator initialization into its module
207!
[2119]208! 2118 2017-01-17 16:38:49Z raasch
209! OpenACC directives removed
210!
[2038]211! 2037 2016-10-26 11:15:40Z knoop
212! Anelastic approximation implemented
213!
[2032]214! 2031 2016-10-21 15:11:58Z knoop
215! renamed variable rho to rho_ocean
216!
[2012]217! 2011 2016-09-19 17:29:57Z kanani
218! Flag urban_surface is now defined in module control_parameters.
219!
[2008]220! 2007 2016-08-24 15:47:17Z kanani
221! Added support for urban surface model,
222! adjusted location_message in case of plant_canopy
223!
[2001]224! 2000 2016-08-20 18:09:15Z knoop
225! Forced header and separation lines into 80 columns
226!
[1993]227! 1992 2016-08-12 15:14:59Z suehring
228! Initializaton of scalarflux at model top
229! Bugfixes in initialization of surface and top salinity flux, top scalar and
230! humidity fluxes
231!
[1961]232! 1960 2016-07-12 16:34:24Z suehring
[1960]233! Separate humidity and passive scalar
234! Increase dimension for mean_inflow_profiles
235! Remove inadvertent write-statement
236! Bugfix, large-scale forcing is still not implemented for passive scalars
[1919]237!
[1958]238! 1957 2016-07-07 10:43:48Z suehring
239! flight module added
240!
[1921]241! 1920 2016-05-30 10:50:15Z suehring
242! Initialize us with very small number to avoid segmentation fault during
243! calculation of Obukhov length
244!
[1919]245! 1918 2016-05-27 14:35:57Z raasch
246! intermediate_timestep_count is set 0 instead 1 for first call of pres,
247! bugfix: initialization of local sum arrays are moved to the beginning of the
248!         routine because otherwise results from pres are overwritten
249!
[1917]250! 1914 2016-05-26 14:44:07Z witha
251! Added initialization of the wind turbine model
252!
[1879]253! 1878 2016-04-19 12:30:36Z hellstea
254! The zeroth element of weight_pres removed as unnecessary
255!
[1851]256! 1849 2016-04-08 11:33:18Z hoffmann
[1849]257! Adapted for modularization of microphysics.
258! precipitation_amount, precipitation_rate, prr moved to arrays_3d.
[1852]259! Initialization of nc_1d, nr_1d, pt_1d, qc_1d, qr_1d, q_1d moved to
[3274]260! bcm_init.
[1849]261!
[1846]262! 1845 2016-04-08 08:29:13Z raasch
263! nzb_2d replaced by nzb_u|v_inner
[1914]264!
[1834]265! 1833 2016-04-07 14:23:03Z raasch
266! initialization of spectra quantities moved to spectra_mod
267!
[1832]268! 1831 2016-04-07 13:15:51Z hoffmann
269! turbulence renamed collision_turbulence
270!
[1827]271! 1826 2016-04-07 12:01:39Z maronga
272! Renamed radiation calls.
273! Renamed canopy model calls.
274!
[1823]275! 1822 2016-04-07 07:49:42Z hoffmann
276! icloud_scheme replaced by microphysics_*
[1914]277!
[1818]278! 1817 2016-04-06 15:44:20Z maronga
279! Renamed lsm calls.
280!
[1816]281! 1815 2016-04-06 13:49:59Z raasch
282! zero-settings for velocities inside topography re-activated (was deactivated
283! in r1762)
284!
[1789]285! 1788 2016-03-10 11:01:04Z maronga
286! Added z0q.
287! Syntax layout improved.
288!
[1784]289! 1783 2016-03-06 18:36:17Z raasch
290! netcdf module name changed + related changes
291!
[1765]292! 1764 2016-02-28 12:45:19Z raasch
293! bugfix: increase size of volume_flow_area_l and volume_flow_initial_l by 1
294!
[1763]295! 1762 2016-02-25 12:31:13Z hellstea
296! Introduction of nested domain feature
297!
[1739]298! 1738 2015-12-18 13:56:05Z raasch
299! calculate mean surface level height for each statistic region
300!
[1735]301! 1734 2015-12-02 12:17:12Z raasch
302! no initial disturbances in case that the disturbance energy limit has been
303! set zero
304!
[1708]305! 1707 2015-11-02 15:24:52Z maronga
306! Bugfix: transfer of Richardson number from 1D model to Obukhov length caused
307! devision by zero in neutral stratification
308!
[1692]309! 1691 2015-10-26 16:17:44Z maronga
310! Call to init_surface_layer added. rif is replaced by ol and zeta.
311!
[1683]312! 1682 2015-10-07 23:56:08Z knoop
313! Code annotations made doxygen readable
314!
[1616]315! 1615 2015-07-08 18:49:19Z suehring
316! Enable turbulent inflow for passive_scalar and humidity
317!
[1586]318! 1585 2015-04-30 07:05:52Z maronga
319! Initialization of radiation code is now done after LSM initializtion
320!
[1576]321! 1575 2015-03-27 09:56:27Z raasch
322! adjustments for psolver-queries
323!
[1552]324! 1551 2015-03-03 14:18:16Z maronga
[1817]325! Allocation of land surface arrays is now done in the subroutine lsm_init_arrays,
[1552]326! which is part of land_surface_model.
327!
[1508]328! 1507 2014-12-10 12:14:18Z suehring
329! Bugfix: set horizontal velocity components to zero inside topography
330!
[1497]331! 1496 2014-12-02 17:25:50Z maronga
332! Added initialization of the land surface and radiation schemes
333!
[1485]334! 1484 2014-10-21 10:53:05Z kanani
[1484]335! Changes due to new module structure of the plant canopy model:
[1508]336! canopy-related initialization (e.g. lad and canopy_heat_flux) moved to new
337! subroutine init_plant_canopy within the module plant_canopy_model_mod,
338! call of subroutine init_plant_canopy added.
[1341]339!
[1432]340! 1431 2014-07-15 14:47:17Z suehring
341! var_d added, in order to normalize spectra.
342!
[1430]343! 1429 2014-07-15 12:53:45Z knoop
344! Ensemble run capability added to parallel random number generator
345!
[1412]346! 1411 2014-05-16 18:01:51Z suehring
347! Initial horizontal velocity profiles were not set to zero at the first vertical
348! grid level in case of non-cyclic lateral boundary conditions.
349!
[1407]350! 1406 2014-05-16 13:47:01Z raasch
351! bugfix: setting of initial velocities at k=1 to zero not in case of a
352! no-slip boundary condition for uv
353!
[1403]354! 1402 2014-05-09 14:25:13Z raasch
355! location messages modified
356!
[1401]357! 1400 2014-05-09 14:03:54Z knoop
358! Parallel random number generator added
359!
[1385]360! 1384 2014-05-02 14:31:06Z raasch
361! location messages added
362!
[1362]363! 1361 2014-04-16 15:17:48Z hoffmann
364! tend_* removed
365! Bugfix: w_subs is not allocated anymore if it is already allocated
366!
[1360]367! 1359 2014-04-11 17:15:14Z hoffmann
368! module lpm_init_mod added to use statements, because lpm_init has become a
369! module
370!
[1354]371! 1353 2014-04-08 15:21:23Z heinze
372! REAL constants provided with KIND-attribute
373!
[1341]374! 1340 2014-03-25 19:45:13Z kanani
375! REAL constants defined as wp-kind
376!
[1323]377! 1322 2014-03-20 16:38:49Z raasch
378! REAL constants defined as wp-kind
379! module interfaces removed
380!
[1321]381! 1320 2014-03-20 08:40:49Z raasch
382! ONLY-attribute added to USE-statements,
383! kind-parameters added to all INTEGER and REAL declaration statements,
384! kinds are defined in new module kinds,
385! revision history before 2012 removed,
386! comment fields (!:) to be used for variable explanations added to
387! all variable declaration statements
388!
[1317]389! 1316 2014-03-17 07:44:59Z heinze
390! Bugfix: allocation of w_subs
391!
[1300]392! 1299 2014-03-06 13:15:21Z heinze
393! Allocate w_subs due to extension of large scale subsidence in combination
394! with large scale forcing data (LSF_DATA)
395!
[1242]396! 1241 2013-10-30 11:36:58Z heinze
397! Overwrite initial profiles in case of nudging
398! Inititialize shf and qsws in case of large_scale_forcing
399!
[1222]400! 1221 2013-09-10 08:59:13Z raasch
401! +rflags_s_inner in copyin statement, use copyin for most arrays instead of
402! copy
403!
[1213]404! 1212 2013-08-15 08:46:27Z raasch
405! array tri is allocated and included in data copy statement
406!
[1196]407! 1195 2013-07-01 12:27:57Z heinze
408! Bugfix: move allocation of ref_state to parin.f90 and read_var_list.f90
409!
[1182]410! 1179 2013-06-14 05:57:58Z raasch
411! allocate and set ref_state to be used in buoyancy terms
412!
[1172]413! 1171 2013-05-30 11:27:45Z raasch
414! diss array is allocated with full size if accelerator boards are used
415!
[1160]416! 1159 2013-05-21 11:58:22Z fricke
417! -bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir
418!
[1154]419! 1153 2013-05-10 14:33:08Z raasch
420! diss array is allocated with dummy elements even if it is not needed
[1171]421! (required by PGI 13.4 / CUDA 5.0)
[1154]422!
[1116]423! 1115 2013-03-26 18:16:16Z hoffmann
424! unused variables removed
425!
[1114]426! 1113 2013-03-10 02:48:14Z raasch
427! openACC directive modified
428!
[1112]429! 1111 2013-03-08 23:54:10Z raasch
430! openACC directives added for pres
431! array diss allocated only if required
432!
[1093]433! 1092 2013-02-02 11:24:22Z raasch
434! unused variables removed
435!
[1066]436! 1065 2012-11-22 17:42:36Z hoffmann
437! allocation of diss (dissipation rate) in case of turbulence = .TRUE. added
438!
[1054]439! 1053 2012-11-13 17:11:03Z hoffmann
[1053]440! allocation and initialisation of necessary data arrays for the two-moment
441! cloud physics scheme the two new prognostic equations (nr, qr):
442! +dr, lambda_r, mu_r, sed_*, xr, *s, *sws, *swst, *, *_p, t*_m, *_1, *_2, *_3,
443! +tend_*, prr
[979]444!
[1037]445! 1036 2012-10-22 13:43:42Z raasch
446! code put under GPL (PALM 3.9)
447!
[1033]448! 1032 2012-10-21 13:03:21Z letzel
449! save memory by not allocating pt_2 in case of neutral = .T.
450!
[1026]451! 1025 2012-10-07 16:04:41Z letzel
452! bugfix: swap indices of mask for ghost boundaries
453!
[1017]454! 1015 2012-09-27 09:23:24Z raasch
455! mask is set to zero for ghost boundaries
456!
[1011]457! 1010 2012-09-20 07:59:54Z raasch
458! cpp switch __nopointer added for pointer free version
459!
[1004]460! 1003 2012-09-14 14:35:53Z raasch
461! nxra,nyna, nzta replaced ny nxr, nyn, nzt
462!
[1002]463! 1001 2012-09-13 14:08:46Z raasch
464! all actions concerning leapfrog scheme removed
465!
[997]466! 996 2012-09-07 10:41:47Z raasch
467! little reformatting
468!
[979]469! 978 2012-08-09 08:28:32Z fricke
[978]470! outflow damping layer removed
471! roughness length for scalar quantites z0h added
472! damping zone for the potential temperatur in case of non-cyclic lateral
473! boundaries added
474! initialization of ptdf_x, ptdf_y
475! initialization of c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l
[708]476!
[850]477! 849 2012-03-15 10:35:09Z raasch
478! init_particles renamed lpm_init
479!
[826]480! 825 2012-02-19 03:03:44Z raasch
481! wang_collision_kernel renamed wang_kernel
482!
[1]483! Revision 1.1  1998/03/09 16:22:22  raasch
484! Initial revision
485!
486!
487! Description:
488! ------------
[1682]489!> Allocation of arrays and initialization of the 3D model via
490!> a) pre-run the 1D model
491!> or
492!> b) pre-set constant linear profiles
493!> or
494!> c) read values of a previous run
[1]495!------------------------------------------------------------------------------!
[1682]496 SUBROUTINE init_3d_model
497 
[1]498
[667]499    USE advec_ws
[1320]500
[1]501    USE arrays_3d
[1849]502
[3274]503    USE basic_constants_and_equations_mod,                                     &
504        ONLY:  c_p, g, l_v, pi, r_d, exner_function, exner_function_invers,    &
505               ideal_gas_law_rho, ideal_gas_law_rho_pt, barometric_formula
506
507    USE bulk_cloud_model_mod,                                                  &
508        ONLY:  bulk_cloud_model, bcm_init, bcm_init_arrays
509
[2696]510    USE chemistry_model_mod,                                                   &
511        ONLY:  chem_emissions
512
[1]513    USE control_parameters
[1320]514   
[1957]515    USE flight_mod,                                                            &
516        ONLY:  flight_init
517   
[1320]518    USE grid_variables,                                                        &
[2037]519        ONLY:  dx, dy, ddx2_mg, ddy2_mg
[2817]520
521    USE gust_mod,                                                              &
522        ONLY:  gust_init, gust_init_arrays, gust_module_enabled
[1320]523   
[1]524    USE indices
[1320]525   
526    USE kinds
[1496]527
528    USE land_surface_model_mod,                                                &
[2232]529        ONLY:  lsm_init, lsm_init_arrays
[3159]530
531    USE lpm_init_mod,                                                          &
532        ONLY:  lpm_init
[1496]533 
[2320]534    USE lsf_nudging_mod,                                                       &
[2696]535        ONLY:  lsf_init, ls_forcing_surf, nudge_init
[1849]536
[2338]537    USE model_1d_mod,                                                          &
[3241]538        ONLY:  init_1d_model, l1d, u1d, v1d
[2338]539
[3159]540    USE multi_agent_system_mod,                                                &
541        ONLY:  agents_active, mas_init
542
[1783]543    USE netcdf_interface,                                                      &
[2817]544        ONLY:  dots_max, dots_num, dots_unit, dots_label
[2696]545
[2906]546    USE netcdf_data_input_mod,                                                 &
[3182]547        ONLY:  init_3d, netcdf_data_input_init_3d
[1320]548   
549    USE particle_attributes,                                                   &
[3241]550        ONLY:  particle_advection
[1320]551   
[1]552    USE pegrid
[1320]553   
[1484]554    USE plant_canopy_model_mod,                                                &
[2746]555        ONLY:  pcm_init
[1496]556
[2934]557    USE pmc_interface,                                                         &
558        ONLY:  nested_run
559
[1496]560    USE radiation_model_mod,                                                   &
[2977]561        ONLY:  average_radiation,                                              &
[2995]562               radiation_init, radiation, radiation_scheme,                    &
[2906]563               radiation_calc_svf, radiation_write_svf,                        &
[2696]564               radiation_interaction, radiation_interactions,                  &
[2920]565               radiation_interaction_init, radiation_read_svf,                 &
[2977]566               radiation_presimulate_solar_pos, radiation_interactions_on
[1484]567   
[1320]568    USE random_function_mod 
569   
[1400]570    USE random_generator_parallel,                                             &
[2172]571        ONLY:  init_parallel_random_generator
[2894]572
573    USE read_restart_data_mod,                                                 &
574        ONLY:  rrd_read_parts_of_global, rrd_local                                     
[1400]575   
[1320]576    USE statistics,                                                            &
[1738]577        ONLY:  hom, hom_sum, mean_surface_level_height, pr_palm, rmask,        &
[1833]578               statistic_regions, sums, sums_divnew_l, sums_divold_l, sums_l,  &
[2277]579               sums_l_l, sums_wsts_bc_l, ts_value,                             &
[1833]580               weight_pres, weight_substep
[2259]581
582    USE synthetic_turbulence_generator_mod,                                    &
[2776]583        ONLY:  stg_init, use_syn_turb_gen
[2259]584
[1691]585    USE surface_layer_fluxes_mod,                                              &
586        ONLY:  init_surface_layer_fluxes
[2232]587
588    USE surface_mod,                                                           &
589        ONLY :  init_surface_arrays, init_surfaces, surf_def_h, surf_lsm_h,    &
[2977]590                surf_usm_h, get_topography_top_index_ji, vertical_surfaces_exist
[1691]591   
[2007]592    USE transpose_indices
[1]593
[2696]594    USE turbulence_closure_mod,                                                &
595        ONLY:  tcm_init_arrays, tcm_init
596
[2007]597    USE urban_surface_mod,                                                     &
[2696]598        ONLY:  usm_init_urban_surface, usm_allocate_surface
[2007]599
[2696]600    USE uv_exposure_model_mod,                                                 &
601        ONLY:  uvem_init, uvem_init_arrays
602
[1914]603    USE wind_turbine_model_mod,                                                &
[2564]604        ONLY:  wtm_init, wtm_init_arrays
[1914]605
[1]606    IMPLICIT NONE
607
[1682]608    INTEGER(iwp) ::  i             !<
609    INTEGER(iwp) ::  ind_array(1)  !<
610    INTEGER(iwp) ::  j             !<
611    INTEGER(iwp) ::  k             !<
[2232]612    INTEGER(iwp) ::  k_surf        !< surface level index
613    INTEGER(iwp) ::  m             !< index of surface element in surface data type
614    INTEGER(iwp) ::  sr            !< index of statistic region
[1]615
[1682]616    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  ngp_2dh_l  !<
[1]617
[1682]618    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l    !<
619    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner_l  !<
[1]620
[2037]621    REAL(wp)     ::  t_surface !< air temperature at the surface
622
[3182]623    REAL(wp), DIMENSION(:), ALLOCATABLE ::  init_l        !< dummy array used for averaging 3D data to obtain inital profiles
[2037]624    REAL(wp), DIMENSION(:), ALLOCATABLE ::  p_hydrostatic !< hydrostatic pressure
625
626    INTEGER(iwp) ::  l       !< loop variable
627    INTEGER(iwp) ::  nzt_l   !< index of top PE boundary for multigrid level
628    REAL(wp) ::  dx_l !< grid spacing along x on different multigrid level
629    REAL(wp) ::  dy_l !< grid spacing along y on different multigrid level
630
[1764]631    REAL(wp), DIMENSION(1:3) ::  volume_flow_area_l     !<
632    REAL(wp), DIMENSION(1:3) ::  volume_flow_initial_l  !<
[1]633
[1738]634    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height_l    !<
[1682]635    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l    !<
636    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_tmp  !<
[1]637
[2550]638    INTEGER(iwp) ::  nz_u_shift   !<
639    INTEGER(iwp) ::  nz_v_shift   !<
640    INTEGER(iwp) ::  nz_w_shift   !<
641    INTEGER(iwp) ::  nz_s_shift   !<
642    INTEGER(iwp) ::  nz_u_shift_l !<
643    INTEGER(iwp) ::  nz_v_shift_l !<
644    INTEGER(iwp) ::  nz_w_shift_l !<
645    INTEGER(iwp) ::  nz_s_shift_l !<
[485]646
[1402]647    CALL location_message( 'allocating arrays', .FALSE. )
[1]648!
649!-- Allocate arrays
[1788]650    ALLOCATE( mean_surface_level_height(0:statistic_regions),                  &
651              mean_surface_level_height_l(0:statistic_regions),                &
652              ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions),    &
653              ngp_3d(0:statistic_regions),                                     &
654              ngp_3d_inner(0:statistic_regions),                               &
655              ngp_3d_inner_l(0:statistic_regions),                             &
656              ngp_3d_inner_tmp(0:statistic_regions),                           &
657              sums_divnew_l(0:statistic_regions),                              &
[1]658              sums_divold_l(0:statistic_regions) )
[1195]659    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
[1788]660    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                    &
661              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),                  &
662              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),                  &
663              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),                &
664              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),                  &
665              sums(nzb:nzt+1,pr_palm+max_pr_user),                             &
666              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),      &
667              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1),    &
668              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                   &
[394]669              ts_value(dots_max,0:statistic_regions) )
[978]670    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
[1]671
[1788]672    ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr),                                    &
673              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
[1010]674              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
675
676#if defined( __nopointer )
[2696]677    ALLOCATE( pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
[1788]678              pt_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
679              u(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
680              u_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
681              v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
682              v_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
683              w(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
684              w_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
685              tpt_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
686              tu_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
687              tv_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
[1010]688              tw_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
689#else
[2696]690    ALLOCATE( pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
[1788]691              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
692              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
693              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
694              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
695              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
696              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
697              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
698              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
699              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
[667]700              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1788]701    IF (  .NOT.  neutral )  THEN
[1032]702       ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
703    ENDIF
[1010]704#endif
705
[673]706!
[707]707!-- Following array is required for perturbation pressure within the iterative
708!-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds
709!-- the weighted average of the substeps and cannot be used in the Poisson
710!-- solver.
711    IF ( psolver == 'sor' )  THEN
712       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1575]713    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
[707]714!
715!--    For performance reasons, multigrid is using one ghost layer only
716       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
[673]717    ENDIF
[1]718
[1111]719!
720!-- Array for storing constant coeffficients of the tridiagonal solver
721    IF ( psolver == 'poisfft' )  THEN
[1212]722       ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) )
[1111]723       ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) )
724    ENDIF
725
[1960]726    IF ( humidity )  THEN
[1]727!
[1960]728!--    3D-humidity
[1010]729#if defined( __nopointer )
[1788]730       ALLOCATE( q(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
731                 q_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
[3011]732                 tq_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
733                 vpt(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1010]734#else
[1788]735       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
736                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
[3011]737                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
738                 vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) 
[1010]739#endif
[1]740
[3011]741       IF ( cloud_droplets )  THEN
[1]742!
[3011]743!--       Liquid water content, change in liquid water content
[1010]744#if defined( __nopointer )
[3011]745          ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
746                     ql_c(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1010]747#else
[3011]748          ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                      &
749                     ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1010]750#endif
751!
[3011]752!--       Real volume of particles (with weighting), volume of particles
753          ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                      &
754                     ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1]755       ENDIF
756
[3011]757    ENDIF   
[1960]758   
759    IF ( passive_scalar )  THEN
[1]760
[1960]761!
762!--    3D scalar arrays
763#if defined( __nopointer )
764       ALLOCATE( s(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
765                 s_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
766                 ts_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
767#else
768       ALLOCATE( s_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
769                 s_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
770                 s_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
771#endif
772    ENDIF
773
[94]774    IF ( ocean )  THEN
[1010]775#if defined( __nopointer )
[1788]776       ALLOCATE( prho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
[2031]777                 rho_ocean(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
[1788]778                 sa(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
779                 sa_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
[1010]780                 tsa_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
781#else
[1788]782       ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
783                 rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                         &
784                 sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
785                 sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
[667]786                 sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[388]787       prho => prho_1
[2031]788       rho_ocean  => rho_1  ! routines calc_mean_profile and diffusion_e require
[388]789                      ! density to be apointer
[1010]790#endif
[94]791    ENDIF
792
[1]793!
[2037]794!-- Allocation of anelastic and Boussinesq approximation specific arrays
795    ALLOCATE( p_hydrostatic(nzb:nzt+1) )
796    ALLOCATE( rho_air(nzb:nzt+1) )
797    ALLOCATE( rho_air_zw(nzb:nzt+1) )
798    ALLOCATE( drho_air(nzb:nzt+1) )
799    ALLOCATE( drho_air_zw(nzb:nzt+1) )
800!
801!-- Density profile calculation for anelastic approximation
[3274]802    t_surface = pt_surface * ( surface_pressure / 1000.0_wp )**( r_d / c_p )
[2037]803    IF ( TRIM( approximation ) == 'anelastic' ) THEN
804       DO  k = nzb, nzt+1
805          p_hydrostatic(k)    = surface_pressure * 100.0_wp *                  &
[3274]806                                ( 1 - ( g * zu(k) ) / ( c_p * t_surface )      &
807                                )**( c_p / r_d )
[2037]808          rho_air(k)          = ( p_hydrostatic(k) *                           &
809                                  ( 100000.0_wp / p_hydrostatic(k)             &
[3274]810                                  )**( r_d / c_p )                             &
[2037]811                                ) / ( r_d * pt_init(k) )
812       ENDDO
813       DO  k = nzb, nzt
814          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
815       ENDDO
816       rho_air_zw(nzt+1)  = rho_air_zw(nzt)                                    &
817                            + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
818    ELSE
[2252]819       DO  k = nzb, nzt+1
820          p_hydrostatic(k)    = surface_pressure * 100.0_wp *                  &
[3274]821                                ( 1 - ( g * zu(nzb) ) / ( c_p * t_surface )    &
822                                )**( c_p / r_d )
[2252]823          rho_air(k)          = ( p_hydrostatic(k) *                           &
824                                  ( 100000.0_wp / p_hydrostatic(k)             &
[3274]825                                  )**( r_d / c_p )                             &
[2252]826                                ) / ( r_d * pt_init(nzb) )
827       ENDDO
828       DO  k = nzb, nzt
829          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
830       ENDDO
831       rho_air_zw(nzt+1)  = rho_air_zw(nzt)                                    &
832                            + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
[2037]833    ENDIF
[2696]834!
[2037]835!-- compute the inverse density array in order to avoid expencive divisions
836    drho_air    = 1.0_wp / rho_air
837    drho_air_zw = 1.0_wp / rho_air_zw
838
839!
840!-- Allocation of flux conversion arrays
841    ALLOCATE( heatflux_input_conversion(nzb:nzt+1) )
842    ALLOCATE( waterflux_input_conversion(nzb:nzt+1) )
843    ALLOCATE( momentumflux_input_conversion(nzb:nzt+1) )
844    ALLOCATE( heatflux_output_conversion(nzb:nzt+1) )
845    ALLOCATE( waterflux_output_conversion(nzb:nzt+1) )
846    ALLOCATE( momentumflux_output_conversion(nzb:nzt+1) )
847
848!
849!-- calculate flux conversion factors according to approximation and in-/output mode
850    DO  k = nzb, nzt+1
851
852        IF ( TRIM( flux_input_mode ) == 'kinematic' )  THEN
853            heatflux_input_conversion(k)      = rho_air_zw(k)
854            waterflux_input_conversion(k)     = rho_air_zw(k)
855            momentumflux_input_conversion(k)  = rho_air_zw(k)
856        ELSEIF ( TRIM( flux_input_mode ) == 'dynamic' ) THEN
[3274]857            heatflux_input_conversion(k)      = 1.0_wp / c_p
[2037]858            waterflux_input_conversion(k)     = 1.0_wp / l_v
859            momentumflux_input_conversion(k)  = 1.0_wp
860        ENDIF
861
862        IF ( TRIM( flux_output_mode ) == 'kinematic' )  THEN
863            heatflux_output_conversion(k)     = drho_air_zw(k)
864            waterflux_output_conversion(k)    = drho_air_zw(k)
865            momentumflux_output_conversion(k) = drho_air_zw(k)
866        ELSEIF ( TRIM( flux_output_mode ) == 'dynamic' ) THEN
[3274]867            heatflux_output_conversion(k)     = c_p
[2037]868            waterflux_output_conversion(k)    = l_v
869            momentumflux_output_conversion(k) = 1.0_wp
870        ENDIF
871
872        IF ( .NOT. humidity ) THEN
873            waterflux_input_conversion(k)  = 1.0_wp
874            waterflux_output_conversion(k) = 1.0_wp
875        ENDIF
876
877    ENDDO
878
879!
880!-- In case of multigrid method, compute grid lengths and grid factors for the
881!-- grid levels with respective density on each grid
882    IF ( psolver(1:9) == 'multigrid' )  THEN
883
884       ALLOCATE( ddx2_mg(maximum_grid_level) )
885       ALLOCATE( ddy2_mg(maximum_grid_level) )
886       ALLOCATE( dzu_mg(nzb+1:nzt+1,maximum_grid_level) )
887       ALLOCATE( dzw_mg(nzb+1:nzt+1,maximum_grid_level) )
888       ALLOCATE( f1_mg(nzb+1:nzt,maximum_grid_level) )
889       ALLOCATE( f2_mg(nzb+1:nzt,maximum_grid_level) )
890       ALLOCATE( f3_mg(nzb+1:nzt,maximum_grid_level) )
891       ALLOCATE( rho_air_mg(nzb:nzt+1,maximum_grid_level) )
892       ALLOCATE( rho_air_zw_mg(nzb:nzt+1,maximum_grid_level) )
893
894       dzu_mg(:,maximum_grid_level) = dzu
895       rho_air_mg(:,maximum_grid_level) = rho_air
896!       
897!--    Next line to ensure an equally spaced grid.
898       dzu_mg(1,maximum_grid_level) = dzu(2)
899       rho_air_mg(nzb,maximum_grid_level) = rho_air(nzb) +                     &
900                                             (rho_air(nzb) - rho_air(nzb+1))
901
902       dzw_mg(:,maximum_grid_level) = dzw
903       rho_air_zw_mg(:,maximum_grid_level) = rho_air_zw
904       nzt_l = nzt
905       DO  l = maximum_grid_level-1, 1, -1
906           dzu_mg(nzb+1,l) = 2.0_wp * dzu_mg(nzb+1,l+1)
907           dzw_mg(nzb+1,l) = 2.0_wp * dzw_mg(nzb+1,l+1)
908           rho_air_mg(nzb,l)    = rho_air_mg(nzb,l+1) + (rho_air_mg(nzb,l+1) - rho_air_mg(nzb+1,l+1))
909           rho_air_zw_mg(nzb,l) = rho_air_zw_mg(nzb,l+1) + (rho_air_zw_mg(nzb,l+1) - rho_air_zw_mg(nzb+1,l+1))
910           rho_air_mg(nzb+1,l)    = rho_air_mg(nzb+1,l+1)
911           rho_air_zw_mg(nzb+1,l) = rho_air_zw_mg(nzb+1,l+1)
912           nzt_l = nzt_l / 2
913           DO  k = 2, nzt_l+1
914              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
915              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
916              rho_air_mg(k,l)    = rho_air_mg(2*k-1,l+1)
917              rho_air_zw_mg(k,l) = rho_air_zw_mg(2*k-1,l+1)
918           ENDDO
919       ENDDO
920
921       nzt_l = nzt
922       dx_l  = dx
923       dy_l  = dy
924       DO  l = maximum_grid_level, 1, -1
925          ddx2_mg(l) = 1.0_wp / dx_l**2
926          ddy2_mg(l) = 1.0_wp / dy_l**2
927          DO  k = nzb+1, nzt_l
928             f2_mg(k,l) = rho_air_zw_mg(k,l) / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
929             f3_mg(k,l) = rho_air_zw_mg(k-1,l) / ( dzu_mg(k,l)   * dzw_mg(k,l) )
930             f1_mg(k,l) = 2.0_wp * ( ddx2_mg(l) + ddy2_mg(l) ) &
931                          * rho_air_mg(k,l) + f2_mg(k,l) + f3_mg(k,l)
932          ENDDO
933          nzt_l = nzt_l / 2
934          dx_l  = dx_l * 2.0_wp
935          dy_l  = dy_l * 2.0_wp
936       ENDDO
937
938    ENDIF
939
940!
[1299]941!-- 1D-array for large scale subsidence velocity
[1361]942    IF ( .NOT. ALLOCATED( w_subs ) )  THEN
943       ALLOCATE ( w_subs(nzb:nzt+1) )
944       w_subs = 0.0_wp
945    ENDIF
[1299]946
947!
[106]948!-- Arrays to store velocity data from t-dt and the phase speeds which
949!-- are needed for radiation boundary conditions
[3182]950    IF ( bc_radiation_l )  THEN
[1788]951       ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2),                               &
952                 v_m_l(nzb:nzt+1,nysg:nyng,0:1),                               &
[667]953                 w_m_l(nzb:nzt+1,nysg:nyng,0:1) )
[73]954    ENDIF
[3182]955    IF ( bc_radiation_r )  THEN
[1788]956       ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
957                 v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
[667]958                 w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) )
[73]959    ENDIF
[3182]960    IF ( bc_radiation_l  .OR.  bc_radiation_r )  THEN
[1788]961       ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng),           &
[667]962                 c_w(nzb:nzt+1,nysg:nyng) )
[106]963    ENDIF
[3182]964    IF ( bc_radiation_s )  THEN
[1788]965       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg),                               &
966                 v_m_s(nzb:nzt+1,1:2,nxlg:nxrg),                               &
[667]967                 w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) )
[73]968    ENDIF
[3182]969    IF ( bc_radiation_n )  THEN
[1788]970       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
971                 v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
[667]972                 w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) )
[73]973    ENDIF
[3182]974    IF ( bc_radiation_s  .OR.  bc_radiation_n )  THEN
[1788]975       ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg),           &
[667]976                 c_w(nzb:nzt+1,nxlg:nxrg) )
[106]977    ENDIF
[3182]978    IF ( bc_radiation_l  .OR.  bc_radiation_r  .OR.  bc_radiation_s  .OR.      &
979         bc_radiation_n )  THEN
[978]980       ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) )                   
981       ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) )
982    ENDIF
[73]983
[978]984
[1010]985#if ! defined( __nopointer )
[73]986!
[1]987!-- Initial assignment of the pointers
[1032]988    IF ( .NOT. neutral )  THEN
989       pt => pt_1;  pt_p => pt_2;  tpt_m => pt_3
990    ELSE
991       pt => pt_1;  pt_p => pt_1;  tpt_m => pt_3
992    ENDIF
[1001]993    u  => u_1;   u_p  => u_2;   tu_m  => u_3
994    v  => v_1;   v_p  => v_2;   tv_m  => v_3
995    w  => w_1;   w_p  => w_2;   tw_m  => w_3
[1]996
[1960]997    IF ( humidity )  THEN
[1001]998       q => q_1;  q_p => q_2;  tq_m => q_3
[3274]999       vpt  => vpt_1
[1001]1000       IF ( cloud_droplets )  THEN
1001          ql   => ql_1
1002          ql_c => ql_2
[1]1003       ENDIF
[1001]1004    ENDIF
[1960]1005   
1006    IF ( passive_scalar )  THEN
1007       s => s_1;  s_p => s_2;  ts_m => s_3
1008    ENDIF   
[1]1009
[1001]1010    IF ( ocean )  THEN
1011       sa => sa_1;  sa_p => sa_2;  tsa_m => sa_3
1012    ENDIF
[1010]1013#endif
[1]1014!
[2696]1015!-- Initialize arrays for turbulence closure
1016    CALL tcm_init_arrays
1017!
1018!-- Initialize surface arrays
[2232]1019    CALL init_surface_arrays
1020!
[3274]1021!-- Allocate microphysics module arrays
1022    IF ( bulk_cloud_model )  THEN
1023       CALL bcm_init_arrays
1024    ENDIF
1025!
[1551]1026!-- Allocate land surface model arrays
1027    IF ( land_surface )  THEN
[1817]1028       CALL lsm_init_arrays
[1551]1029    ENDIF
1030
1031!
[1914]1032!-- Allocate wind turbine model arrays
1033    IF ( wind_turbine )  THEN
1034       CALL wtm_init_arrays
1035    ENDIF
[1957]1036!
[2817]1037!-- Allocate gust module arrays
1038    IF ( gust_module_enabled )  THEN
1039       CALL gust_init_arrays
1040    ENDIF
1041
1042!
[1957]1043!-- Initialize virtual flight measurements
1044    IF ( virtual_flight )  THEN
1045       CALL flight_init
1046    ENDIF
[1914]1047
1048!
[2696]1049!-- Read uv exposure input data
1050    IF ( uv_exposure )  THEN
1051       CALL uvem_init
[2320]1052    ENDIF
1053!
[2696]1054!-- Allocate uv exposure arrays
1055    IF ( uv_exposure )  THEN
1056       CALL uvem_init_arrays
[2320]1057    ENDIF
1058
1059!
[709]1060!-- Allocate arrays containing the RK coefficient for calculation of
1061!-- perturbation pressure and turbulent fluxes. At this point values are
1062!-- set for pressure calculation during initialization (where no timestep
1063!-- is done). Further below the values needed within the timestep scheme
1064!-- will be set.
[1788]1065    ALLOCATE( weight_substep(1:intermediate_timestep_count_max),               &
[1878]1066              weight_pres(1:intermediate_timestep_count_max) )
[1340]1067    weight_substep = 1.0_wp
1068    weight_pres    = 1.0_wp
[1918]1069    intermediate_timestep_count = 0  ! needed when simulated_time = 0.0
[673]1070       
[1402]1071    CALL location_message( 'finished', .TRUE. )
[1918]1072
[673]1073!
[3014]1074!-- Initialize time series
1075    ts_value = 0.0_wp
1076
1077!
[1918]1078!-- Initialize local summation arrays for routine flow_statistics.
1079!-- This is necessary because they may not yet have been initialized when they
1080!-- are called from flow_statistics (or - depending on the chosen model run -
1081!-- are never initialized)
1082    sums_divnew_l      = 0.0_wp
1083    sums_divold_l      = 0.0_wp
1084    sums_l_l           = 0.0_wp
1085    sums_wsts_bc_l     = 0.0_wp
[3182]1086   
[1918]1087!
[1]1088!-- Initialize model variables
[1788]1089    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
[328]1090         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
[1]1091!
[2696]1092!--    Initialization with provided input data derived from larger-scale model
1093       IF ( INDEX( initializing_actions, 'inifor' ) /= 0 )  THEN
1094          CALL location_message( 'initializing with INIFOR', .FALSE. )
1095!
[3051]1096!--       Read initial 1D profiles or 3D data from NetCDF file, depending
1097!--       on the provided level-of-detail.
[2696]1098!--       At the moment, only u, v, w, pt and q are provided.
1099          CALL netcdf_data_input_init_3d
1100!
[3182]1101!--       Please note, Inifor provides data from nzb+1 to nzt.
1102!--       Bottom and top boundary conditions for Inifor profiles are already
1103!--       set (just after reading), so that this is not necessary here.
1104!--       Depending on the provided level-of-detail, initial Inifor data is
1105!--       either stored on data type (lod=1), or directly on 3D arrays (lod=2).
1106!--       In order to obtain also initial profiles in case of lod=2 (which
1107!--       is required for e.g. damping), average over 3D data.
1108          IF( init_3d%lod_u == 1 )  THEN
1109             u_init = init_3d%u_init
1110          ELSEIF( init_3d%lod_u == 2 )  THEN
1111             ALLOCATE( init_l(nzb:nzt+1) ) 
1112             DO  k = nzb, nzt+1
1113                init_l(k) = SUM( u(k,nys:nyn,nxl:nxr) )
1114             ENDDO
1115             init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
[1384]1116
[3182]1117#if defined( __parallel )
1118             CALL MPI_ALLREDUCE( init_l, u_init, nzt+1-nzb+1,                  &
1119                                 MPI_REAL, MPI_SUM, comm2d, ierr )
1120#else
1121             u_init = init_l
1122#endif
1123             DEALLOCATE( init_l )
[3051]1124
[2696]1125          ENDIF
[3182]1126           
1127          IF( init_3d%lod_v == 1 )  THEN 
1128             v_init = init_3d%v_init
1129          ELSEIF( init_3d%lod_v == 2 )  THEN
1130             ALLOCATE( init_l(nzb:nzt+1) ) 
1131             DO  k = nzb, nzt+1
1132                init_l(k) = SUM( v(k,nys:nyn,nxl:nxr) )
1133             ENDDO
1134             init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
[2696]1135
[3182]1136#if defined( __parallel )
1137             CALL MPI_ALLREDUCE( init_l, v_init, nzt+1-nzb+1,                  &
1138                                 MPI_REAL, MPI_SUM, comm2d, ierr )
1139#else
1140             v_init = init_l
1141#endif
1142             DEALLOCATE( init_l )
1143          ENDIF
1144          IF( .NOT. neutral )  THEN
1145             IF( init_3d%lod_pt == 1 )  THEN
1146                pt_init = init_3d%pt_init
1147             ELSEIF( init_3d%lod_pt == 2 )  THEN
1148                ALLOCATE( init_l(nzb:nzt+1) ) 
1149                DO  k = nzb, nzt+1
1150                   init_l(k) = SUM( pt(k,nys:nyn,nxl:nxr) )
1151                ENDDO
1152                init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
1153
1154#if defined( __parallel )
1155                CALL MPI_ALLREDUCE( init_l, pt_init, nzt+1-nzb+1,               &
1156                                    MPI_REAL, MPI_SUM, comm2d, ierr )
1157#else
1158                pt_init = init_l
1159#endif
1160                DEALLOCATE( init_l )
1161             ENDIF
1162          ENDIF
1163
1164
1165          IF( humidity )  THEN
1166             IF( init_3d%lod_q == 1 )  THEN
1167                q_init = init_3d%q_init
1168             ELSEIF( init_3d%lod_q == 2 )  THEN
1169                ALLOCATE( init_l(nzb:nzt+1) ) 
1170                DO  k = nzb, nzt+1
1171                   init_l(k) = SUM( q(k,nys:nyn,nxl:nxr) )
1172                ENDDO
1173                init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
1174
1175#if defined( __parallel )
1176                CALL MPI_ALLREDUCE( init_l, q_init, nzt+1-nzb+1,               &
1177                                    MPI_REAL, MPI_SUM, comm2d, ierr )
1178#else
1179                q_init = init_l
1180#endif
1181                DEALLOCATE( init_l )
1182             ENDIF
1183          ENDIF
1184
[2696]1185!
[3182]1186!--       Write initial profiles onto 3D arrays. Note, only in case of lod = 1,
1187!--       for lod = 2 data is already on 3D arrays.   
[2696]1188          DO  i = nxlg, nxrg
1189             DO  j = nysg, nyng
[3051]1190                IF( init_3d%lod_u == 1 )  u(:,j,i) = u_init(:)
1191                IF( init_3d%lod_v == 1 )  v(:,j,i) = v_init(:)
1192                IF( .NOT. neutral  .AND.  init_3d%lod_pt == 1 )                &
1193                   pt(:,j,i) = pt_init(:)
[3182]1194                IF( humidity  .AND.  init_3d%lod_q == 1 )  q(:,j,i) = q_init(:)
[2696]1195             ENDDO
1196          ENDDO
1197!
[3182]1198!--       Exchange ghost points in case of level-of-detail = 2
1199          IF( init_3d%lod_u == 2 )   CALL exchange_horiz( u, nbgp )
1200          IF( init_3d%lod_v == 2 )   CALL exchange_horiz( v, nbgp )
1201          IF( init_3d%lod_w == 2 )   CALL exchange_horiz( w, nbgp )
1202          IF( .NOT. neutral  .AND.  init_3d%lod_pt == 2 )                      &
1203             CALL exchange_horiz( pt, nbgp )
1204          IF( humidity  .AND.  init_3d%lod_q == 2 )                            &
1205             CALL exchange_horiz( q, nbgp )
1206!
1207!--       Set geostrophic wind components. 
[2938]1208          IF ( init_3d%from_file_ug )  THEN
1209             ug(:) = init_3d%ug_init(:)
1210          ENDIF
1211          IF ( init_3d%from_file_vg )  THEN
1212             vg(:) = init_3d%vg_init(:)
1213          ENDIF
[3182]1214         
[2938]1215          ug(nzt+1) = ug(nzt)
1216          vg(nzt+1) = vg(nzt)
[2696]1217!
1218!--       Set inital w to 0
1219          w = 0.0_wp
1220
1221          IF ( passive_scalar )  THEN
1222             DO  i = nxlg, nxrg
1223                DO  j = nysg, nyng
1224                   s(:,j,i) = s_init
1225                ENDDO
1226             ENDDO
1227          ENDIF
1228
1229          IF ( ocean )  THEN
1230             DO  i = nxlg, nxrg
1231                DO  j = nysg, nyng
1232                   sa(:,j,i) = sa_init
1233                ENDDO
1234             ENDDO
1235          ENDIF
1236!
1237!--       Set velocity components at non-atmospheric / oceanic grid points to
1238!--       zero.
1239          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1240          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )
1241          w = MERGE( w, 0.0_wp, BTEST( wall_flags_0, 3 ) )
[2700]1242!
1243!--       Initialize surface variables, e.g. friction velocity, momentum
1244!--       fluxes, etc.
1245          CALL init_surfaces
[2938]1246!
1247!--       Initialize turbulence generator
1248          IF( use_syn_turb_gen )  CALL stg_init
[2696]1249
1250          CALL location_message( 'finished', .TRUE. )
1251!
1252!--    Initialization via computed 1D-model profiles
1253       ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1254
[1402]1255          CALL location_message( 'initializing with 1D model profiles', .FALSE. )
[1]1256!
1257!--       Use solutions of the 1D model as initial profiles,
1258!--       start 1D model
1259          CALL init_1d_model
1260!
1261!--       Transfer initial profiles to the arrays of the 3D model
[667]1262          DO  i = nxlg, nxrg
1263             DO  j = nysg, nyng
[1]1264                pt(:,j,i) = pt_init
1265                u(:,j,i)  = u1d
1266                v(:,j,i)  = v1d
1267             ENDDO
1268          ENDDO
1269
[1960]1270          IF ( humidity )  THEN
[667]1271             DO  i = nxlg, nxrg
1272                DO  j = nysg, nyng
[1]1273                   q(:,j,i) = q_init
1274                ENDDO
1275             ENDDO
1276          ENDIF
[2292]1277
[1960]1278          IF ( passive_scalar )  THEN
1279             DO  i = nxlg, nxrg
1280                DO  j = nysg, nyng
1281                   s(:,j,i) = s_init
1282                ENDDO
1283             ENDDO   
1284          ENDIF
[1]1285!
1286!--          Store initial profiles for output purposes etc.
[2696]1287          IF ( .NOT. constant_diffusion )  THEN
[1]1288             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
1289          ENDIF
1290!
[2696]1291!--       Set velocities back to zero
[2758]1292          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1293          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )         
[1]1294!
[2696]1295!--       WARNING: The extra boundary conditions set after running the
1296!--       -------  1D model impose an error on the divergence one layer
1297!--                below the topography; need to correct later
1298!--       ATTENTION: Provisional correction for Piacsek & Williams
1299!--       ---------  advection scheme: keep u and v zero one layer below
1300!--                  the topography.
1301          IF ( ibc_uv_b == 1 )  THEN
[667]1302!
[2696]1303!--          Neumann condition
1304             DO  i = nxl-1, nxr+1
1305                DO  j = nys-1, nyn+1
1306                   u(nzb,j,i) = u(nzb+1,j,i)
1307                   v(nzb,j,i) = v(nzb+1,j,i)
[1]1308                ENDDO
[2696]1309             ENDDO
[1]1310
1311          ENDIF
[2618]1312!
1313!--       Initialize surface variables, e.g. friction velocity, momentum
1314!--       fluxes, etc.
1315          CALL init_surfaces
[1]1316
[1402]1317          CALL location_message( 'finished', .TRUE. )
[1384]1318
[1788]1319       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )    &
[1]1320       THEN
[1241]1321
[1402]1322          CALL location_message( 'initializing with constant profiles', .FALSE. )
[1]1323!
[2259]1324!--       Overwrite initial profiles in case of synthetic turbulence generator
[2938]1325          IF( use_syn_turb_gen )  CALL stg_init
[2259]1326
1327!
[1]1328!--       Use constructed initial profiles (velocity constant with height,
1329!--       temperature profile with constant gradient)
[667]1330          DO  i = nxlg, nxrg
1331             DO  j = nysg, nyng
[1]1332                pt(:,j,i) = pt_init
1333                u(:,j,i)  = u_init
1334                v(:,j,i)  = v_init
1335             ENDDO
1336          ENDDO
1337!
[2758]1338!--       Mask topography
1339          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1340          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )
1341!
[292]1342!--       Set initial horizontal velocities at the lowest computational grid
1343!--       levels to zero in order to avoid too small time steps caused by the
1344!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
[2758]1345!--       in the limiting formula!).
1346!--       Please note, in case land- or urban-surface model is used and a
1347!--       spinup is applied, masking the lowest computational level is not
1348!--       possible as MOST as well as energy-balance parametrizations will not
1349!--       work with zero wind velocity.
1350          IF ( ibc_uv_b /= 1  .AND.  .NOT.  spinup )  THEN
[1815]1351             DO  i = nxlg, nxrg
1352                DO  j = nysg, nyng
[2232]1353                   DO  k = nzb, nzt
1354                      u(k,j,i) = MERGE( u(k,j,i), 0.0_wp,                      &
1355                                        BTEST( wall_flags_0(k,j,i), 20 ) )
1356                      v(k,j,i) = MERGE( v(k,j,i), 0.0_wp,                      &
1357                                        BTEST( wall_flags_0(k,j,i), 21 ) )
1358                   ENDDO
[1815]1359                ENDDO
1360             ENDDO
1361          ENDIF
[1]1362
[1960]1363          IF ( humidity )  THEN
[667]1364             DO  i = nxlg, nxrg
1365                DO  j = nysg, nyng
[1]1366                   q(:,j,i) = q_init
1367                ENDDO
1368             ENDDO
1369          ENDIF
[1960]1370         
1371          IF ( passive_scalar )  THEN
1372             DO  i = nxlg, nxrg
1373                DO  j = nysg, nyng
1374                   s(:,j,i) = s_init
1375                ENDDO
1376             ENDDO
1377          ENDIF
[1]1378
[94]1379          IF ( ocean )  THEN
[667]1380             DO  i = nxlg, nxrg
1381                DO  j = nysg, nyng
[94]1382                   sa(:,j,i) = sa_init
1383                ENDDO
1384             ENDDO
1385          ENDIF
[1920]1386!
[1]1387!--       Compute initial temperature field and other constants used in case
1388!--       of a sloping surface
1389          IF ( sloping_surface )  CALL init_slope
[2618]1390!
1391!--       Initialize surface variables, e.g. friction velocity, momentum
1392!--       fluxes, etc.
1393          CALL init_surfaces
[1]1394
[1402]1395          CALL location_message( 'finished', .TRUE. )
[1384]1396
[1788]1397       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )                  &
[46]1398       THEN
[1384]1399
[1402]1400          CALL location_message( 'initializing by user', .FALSE. )
[46]1401!
[2618]1402!--       Pre-initialize surface variables, i.e. setting start- and end-indices
1403!--       at each (j,i)-location. Please note, this does not supersede
1404!--       user-defined initialization of surface quantities.
1405          CALL init_surfaces
1406!
[46]1407!--       Initialization will completely be done by the user
1408          CALL user_init_3d_model
1409
[1402]1410          CALL location_message( 'finished', .TRUE. )
[1384]1411
[1]1412       ENDIF
[1384]1413
[1402]1414       CALL location_message( 'initializing statistics, boundary conditions, etc.', &
1415                              .FALSE. )
[1384]1416
[667]1417!
1418!--    Bottom boundary
1419       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
[1340]1420          u(nzb,:,:) = 0.0_wp
1421          v(nzb,:,:) = 0.0_wp
[667]1422       ENDIF
[1]1423
1424!
[151]1425!--    Apply channel flow boundary condition
[132]1426       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
[1340]1427          u(nzt+1,:,:) = 0.0_wp
1428          v(nzt+1,:,:) = 0.0_wp
[132]1429       ENDIF
1430
1431!
[1]1432!--    Calculate virtual potential temperature
[1960]1433       IF ( humidity )  vpt = pt * ( 1.0_wp + 0.61_wp * q )
[1]1434
1435!
[2696]1436!--    Store initial profiles for output purposes etc.. Please note, in case of
1437!--    initialization of u, v, w, pt, and q via output data derived from larger
1438!--    scale models, data will not be horizontally homogeneous. Actually, a mean
1439!--    profile should be calculated before.   
[1]1440       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
1441       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
[667]1442       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
[1340]1443          hom(nzb,1,5,:) = 0.0_wp
1444          hom(nzb,1,6,:) = 0.0_wp
[1]1445       ENDIF
1446       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1447
[2696]1448
1449!
1450!--    Store initial salinity profile
[97]1451       IF ( ocean )  THEN
1452          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
1453       ENDIF
[1]1454
[75]1455       IF ( humidity )  THEN
[1]1456!
1457!--       Store initial profile of total water content, virtual potential
1458!--       temperature
1459          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
1460          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
[2696]1461!
[3040]1462!--       Store initial profile of mixing ratio and potential
[2696]1463!--       temperature
[3274]1464          IF ( bulk_cloud_model  .OR.  cloud_droplets ) THEN
[1]1465             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1466             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1467          ENDIF
1468       ENDIF
1469
[2696]1470!
1471!--    Store initial scalar profile
[1]1472       IF ( passive_scalar )  THEN
[2513]1473          hom(:,1,121,:) = SPREAD(  s(:,nys,nxl), 2, statistic_regions+1 )
[1]1474       ENDIF
1475
1476!
[1400]1477!--    Initialize the random number generators (from numerical recipes)
1478       CALL random_function_ini
[1429]1479       
[1400]1480       IF ( random_generator == 'random-parallel' )  THEN
[3241]1481          CALL init_parallel_random_generator( nx, nys, nyn, nxl, nxr )
[1400]1482       ENDIF
1483!
[1179]1484!--    Set the reference state to be used in the buoyancy terms (for ocean runs
1485!--    the reference state will be set (overwritten) in init_ocean)
1486       IF ( use_single_reference_value )  THEN
[1788]1487          IF (  .NOT.  humidity )  THEN
[1179]1488             ref_state(:) = pt_reference
1489          ELSE
1490             ref_state(:) = vpt_reference
1491          ENDIF
1492       ELSE
[1788]1493          IF (  .NOT.  humidity )  THEN
[1179]1494             ref_state(:) = pt_init(:)
1495          ELSE
1496             ref_state(:) = vpt(:,nys,nxl)
1497          ENDIF
1498       ENDIF
[152]1499
1500!
[707]1501!--    For the moment, vertical velocity is zero
[1340]1502       w = 0.0_wp
[1]1503
1504!
1505!--    Initialize array sums (must be defined in first call of pres)
[1340]1506       sums = 0.0_wp
[1]1507
1508!
[707]1509!--    In case of iterative solvers, p must get an initial value
[1575]1510       IF ( psolver(1:9) == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0_wp
[707]1511!
[1]1512!--    Impose vortex with vertical axis on the initial velocity profile
1513       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
1514          CALL init_rankine
1515       ENDIF
1516
1517!
[3035]1518!--    Impose temperature anomaly (advection test only) or warm air bubble
1519!--    close to surface
1520       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0  .OR.  &
1521            INDEX( initializing_actions, 'initialize_bubble' ) /= 0  )  THEN
[1]1522          CALL init_pt_anomaly
1523       ENDIF
[3035]1524       
[1]1525!
1526!--    If required, change the surface temperature at the start of the 3D run
[1340]1527       IF ( pt_surface_initial_change /= 0.0_wp )  THEN
[1]1528          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
1529       ENDIF
1530
1531!
1532!--    If required, change the surface humidity/scalar at the start of the 3D
1533!--    run
[1960]1534       IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )              &
[1]1535          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
[1960]1536         
1537       IF ( passive_scalar .AND.  s_surface_initial_change /= 0.0_wp )         &
1538          s(nzb,:,:) = s(nzb,:,:) + s_surface_initial_change
1539       
[1]1540
1541!
1542!--    Initialize old and new time levels.
[2696]1543       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1544       pt_p = pt; u_p = u; v_p = v; w_p = w
[1]1545
[1960]1546       IF ( humidity  )  THEN
[1340]1547          tq_m = 0.0_wp
[1]1548          q_p = q
1549       ENDIF
[1960]1550       
1551       IF ( passive_scalar )  THEN
1552          ts_m = 0.0_wp
1553          s_p  = s
1554       ENDIF       
[1]1555
[94]1556       IF ( ocean )  THEN
[1340]1557          tsa_m = 0.0_wp
[94]1558          sa_p  = sa
1559       ENDIF
[667]1560       
[1402]1561       CALL location_message( 'finished', .TRUE. )
[94]1562
[1788]1563    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.         &
[2232]1564             TRIM( initializing_actions ) == 'cyclic_fill' )                   &
[1]1565    THEN
[1384]1566
[1402]1567       CALL location_message( 'initializing in case of restart / cyclic_fill', &
1568                              .FALSE. )
[1]1569!
[2232]1570!--    Initialize surface elements and its attributes, e.g. heat- and
1571!--    momentumfluxes, roughness, scaling parameters. As number of surface
1572!--    elements might be different between runs, e.g. in case of cyclic fill,
1573!--    and not all surface elements are read, surface elements need to be
1574!--    initialized before.     
1575       CALL init_surfaces
1576!
[767]1577!--    When reading data for cyclic fill of 3D prerun data files, read
1578!--    some of the global variables from the restart file which are required
1579!--    for initializing the inflow
[328]1580       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[559]1581
[759]1582          DO  i = 0, io_blocks-1
1583             IF ( i == io_group )  THEN
[2894]1584                CALL rrd_read_parts_of_global
[759]1585             ENDIF
1586#if defined( __parallel )
1587             CALL MPI_BARRIER( comm2d, ierr )
1588#endif
1589          ENDDO
[328]1590
[767]1591       ENDIF
1592
[151]1593!
[2894]1594!--    Read processor specific binary data from restart file
[767]1595       DO  i = 0, io_blocks-1
1596          IF ( i == io_group )  THEN
[2894]1597             CALL rrd_local
[767]1598          ENDIF
1599#if defined( __parallel )
1600          CALL MPI_BARRIER( comm2d, ierr )
1601#endif
1602       ENDDO
1603
[328]1604!
[2550]1605!--    In case of complex terrain and cyclic fill method as initialization,
1606!--    shift initial data in the vertical direction for each point in the
1607!--    x-y-plane depending on local surface height
1608       IF ( complex_terrain  .AND.                                             &
1609            TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1610          DO  i = nxlg, nxrg
1611             DO  j = nysg, nyng
[2698]1612                nz_u_shift = get_topography_top_index_ji( j, i, 'u' )
1613                nz_v_shift = get_topography_top_index_ji( j, i, 'v' )
1614                nz_w_shift = get_topography_top_index_ji( j, i, 'w' )
1615                nz_s_shift = get_topography_top_index_ji( j, i, 's' )
[2550]1616
1617                u(nz_u_shift:nzt+1,j,i)  = u(0:nzt+1-nz_u_shift,j,i)               
1618
1619                v(nz_v_shift:nzt+1,j,i)  = v(0:nzt+1-nz_v_shift,j,i)
1620
1621                w(nz_w_shift:nzt+1,j,i)  = w(0:nzt+1-nz_w_shift,j,i)
1622
1623                p(nz_s_shift:nzt+1,j,i)  =  p(0:nzt+1-nz_s_shift,j,i)
1624                pt(nz_s_shift:nzt+1,j,i) = pt(0:nzt+1-nz_s_shift,j,i)
1625             ENDDO
1626          ENDDO
1627       ENDIF
1628
1629!
[767]1630!--    Initialization of the turbulence recycling method
[1788]1631       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.               &
[767]1632            turbulent_inflow )  THEN
1633!
1634!--       First store the profiles to be used at the inflow.
1635!--       These profiles are the (temporally) and horizontally averaged vertical
1636!--       profiles from the prerun. Alternatively, prescribed profiles
1637!--       for u,v-components can be used.
[3288]1638          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,1:num_mean_inflow_profiles) )
[151]1639
[767]1640          IF ( use_prescribed_profile_data )  THEN
1641             mean_inflow_profiles(:,1) = u_init            ! u
1642             mean_inflow_profiles(:,2) = v_init            ! v
1643          ELSE
[328]1644             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
1645             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
[767]1646          ENDIF
1647          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
[1960]1648          IF ( humidity )                                                      &
1649             mean_inflow_profiles(:,6) = hom_sum(:,41,0)   ! q
1650          IF ( passive_scalar )                                                &
1651             mean_inflow_profiles(:,7) = hom_sum(:,115,0)   ! s
[2550]1652!
1653!--       In case of complex terrain, determine vertical displacement at inflow
1654!--       boundary and adjust mean inflow profiles
1655          IF ( complex_terrain )  THEN
1656             IF ( nxlg <= 0 .AND. nxrg >= 0 .AND. nysg <= 0 .AND. nyng >= 0 )  THEN
[2698]1657                nz_u_shift_l = get_topography_top_index_ji( 0, 0, 'u' )
1658                nz_v_shift_l = get_topography_top_index_ji( 0, 0, 'v' )
1659                nz_w_shift_l = get_topography_top_index_ji( 0, 0, 'w' )
1660                nz_s_shift_l = get_topography_top_index_ji( 0, 0, 's' )
[2550]1661             ELSE
1662                nz_u_shift_l = 0
1663                nz_v_shift_l = 0
1664                nz_w_shift_l = 0
1665                nz_s_shift_l = 0
1666             ENDIF
[151]1667
[2550]1668#if defined( __parallel )
1669             CALL MPI_ALLREDUCE(nz_u_shift_l, nz_u_shift, 1, MPI_INTEGER,      &
1670                                MPI_MAX, comm2d, ierr)
1671             CALL MPI_ALLREDUCE(nz_v_shift_l, nz_v_shift, 1, MPI_INTEGER,      &
1672                                MPI_MAX, comm2d, ierr)
1673             CALL MPI_ALLREDUCE(nz_w_shift_l, nz_w_shift, 1, MPI_INTEGER,      & 
1674                                MPI_MAX, comm2d, ierr)
1675             CALL MPI_ALLREDUCE(nz_s_shift_l, nz_s_shift, 1, MPI_INTEGER,      &
1676                                MPI_MAX, comm2d, ierr)
1677#else
1678             nz_u_shift = nz_u_shift_l
1679             nz_v_shift = nz_v_shift_l
1680             nz_w_shift = nz_w_shift_l
1681             nz_s_shift = nz_s_shift_l
1682#endif
1683
1684             mean_inflow_profiles(:,1) = 0.0_wp
1685             mean_inflow_profiles(nz_u_shift:nzt+1,1) = hom_sum(0:nzt+1-nz_u_shift,1,0)  ! u
1686
1687             mean_inflow_profiles(:,2) = 0.0_wp
1688             mean_inflow_profiles(nz_v_shift:nzt+1,2) = hom_sum(0:nzt+1-nz_v_shift,2,0)  ! v
1689
1690             mean_inflow_profiles(nz_s_shift:nzt+1,4) = hom_sum(0:nzt+1-nz_s_shift,4,0)  ! pt
1691
1692          ENDIF
1693
[151]1694!
[767]1695!--       If necessary, adjust the horizontal flow field to the prescribed
1696!--       profiles
1697          IF ( use_prescribed_profile_data )  THEN
1698             DO  i = nxlg, nxrg
[667]1699                DO  j = nysg, nyng
[328]1700                   DO  k = nzb, nzt+1
[767]1701                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1702                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
[328]1703                   ENDDO
[151]1704                ENDDO
[767]1705             ENDDO
1706          ENDIF
[151]1707
1708!
[767]1709!--       Use these mean profiles at the inflow (provided that Dirichlet
1710!--       conditions are used)
[3182]1711          IF ( bc_dirichlet_l )  THEN
[767]1712             DO  j = nysg, nyng
1713                DO  k = nzb, nzt+1
1714                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1715                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
[1340]1716                   w(k,j,nxlg:-1)  = 0.0_wp
[767]1717                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
[1960]1718                   IF ( humidity )                                             &
[1615]1719                      q(k,j,nxlg:-1)  = mean_inflow_profiles(k,6)
[1960]1720                   IF ( passive_scalar )                                       &
1721                      s(k,j,nxlg:-1)  = mean_inflow_profiles(k,7)                     
[767]1722                ENDDO
1723             ENDDO
1724          ENDIF
1725
[151]1726!
[767]1727!--       Calculate the damping factors to be used at the inflow. For a
1728!--       turbulent inflow the turbulent fluctuations have to be limited
1729!--       vertically because otherwise the turbulent inflow layer will grow
1730!--       in time.
[1340]1731          IF ( inflow_damping_height == 9999999.9_wp )  THEN
[767]1732!
1733!--          Default: use the inversion height calculated by the prerun; if
1734!--          this is zero, inflow_damping_height must be explicitly
1735!--          specified.
[1340]1736             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0_wp )  THEN
[767]1737                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1738             ELSE
[1788]1739                WRITE( message_string, * ) 'inflow_damping_height must be ',   &
1740                     'explicitly specified because&the inversion height ',     &
[767]1741                     'calculated by the prerun is zero.'
1742                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
[292]1743             ENDIF
[151]1744
[767]1745          ENDIF
1746
[1340]1747          IF ( inflow_damping_width == 9999999.9_wp )  THEN
[151]1748!
[767]1749!--          Default for the transition range: one tenth of the undamped
1750!--          layer
[1340]1751             inflow_damping_width = 0.1_wp * inflow_damping_height
[151]1752
[767]1753          ENDIF
[151]1754
[767]1755          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
[151]1756
[767]1757          DO  k = nzb, nzt+1
[151]1758
[767]1759             IF ( zu(k) <= inflow_damping_height )  THEN
[1340]1760                inflow_damping_factor(k) = 1.0_wp
[996]1761             ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) )  THEN
[1340]1762                inflow_damping_factor(k) = 1.0_wp -                            &
[996]1763                                           ( zu(k) - inflow_damping_height ) / &
1764                                           inflow_damping_width
[767]1765             ELSE
[1340]1766                inflow_damping_factor(k) = 0.0_wp
[767]1767             ENDIF
[151]1768
[767]1769          ENDDO
[151]1770
[147]1771       ENDIF
1772
[152]1773!
[2696]1774!--    Inside buildings set velocities back to zero
[1788]1775       IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND.                &
[359]1776            topography /= 'flat' )  THEN
1777!
[2696]1778!--       Inside buildings set velocities back to zero.
1779!--       Other scalars (pt, q, s, p, sa, ...) are ignored at present,
[359]1780!--       maybe revise later.
[1001]1781          DO  i = nxlg, nxrg
1782             DO  j = nysg, nyng
[2232]1783                DO  k = nzb, nzt
1784                   u(k,j,i)     = MERGE( u(k,j,i), 0.0_wp,                     &
1785                                         BTEST( wall_flags_0(k,j,i), 1 ) )
1786                   v(k,j,i)     = MERGE( v(k,j,i), 0.0_wp,                     &
1787                                         BTEST( wall_flags_0(k,j,i), 2 ) )
1788                   w(k,j,i)     = MERGE( w(k,j,i), 0.0_wp,                     &
1789                                         BTEST( wall_flags_0(k,j,i), 3 ) )
1790                ENDDO
[359]1791             ENDDO
[1001]1792          ENDDO
[359]1793
1794       ENDIF
1795
1796!
[1]1797!--    Calculate initial temperature field and other constants used in case
1798!--    of a sloping surface
1799       IF ( sloping_surface )  CALL init_slope
1800
1801!
1802!--    Initialize new time levels (only done in order to set boundary values
1803!--    including ghost points)
[2696]1804       pt_p = pt; u_p = u; v_p = v; w_p = w
[1960]1805       IF ( humidity )  THEN
[1053]1806          q_p = q
1807       ENDIF
[1960]1808       IF ( passive_scalar )  s_p  = s
1809       IF ( ocean          )  sa_p = sa
[1]1810
[181]1811!
1812!--    Allthough tendency arrays are set in prognostic_equations, they have
1813!--    have to be predefined here because they are used (but multiplied with 0)
1814!--    there before they are set.
[2696]1815       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
[1960]1816       IF ( humidity )  THEN
[1340]1817          tq_m = 0.0_wp
[1053]1818       ENDIF
[1960]1819       IF ( passive_scalar )  ts_m  = 0.0_wp
1820       IF ( ocean          )  tsa_m = 0.0_wp
[2259]1821!
1822!--    Initialize synthetic turbulence generator in case of restart.
1823       IF ( TRIM( initializing_actions ) == 'read_restart_data'  .AND.         &
[2776]1824            use_syn_turb_gen )  CALL stg_init
[181]1825
[1402]1826       CALL location_message( 'finished', .TRUE. )
[1384]1827
[1]1828    ELSE
1829!
1830!--    Actually this part of the programm should not be reached
[254]1831       message_string = 'unknown initializing problem'
1832       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
[1]1833    ENDIF
1834
[2696]1835!
1836!-- Initialize TKE, Kh and Km
1837    CALL tcm_init
[151]1838
[2696]1839
[151]1840    IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
[1]1841!
[151]1842!--    Initialize old timelevels needed for radiation boundary conditions
[3182]1843       IF ( bc_radiation_l )  THEN
[151]1844          u_m_l(:,:,:) = u(:,:,1:2)
1845          v_m_l(:,:,:) = v(:,:,0:1)
1846          w_m_l(:,:,:) = w(:,:,0:1)
1847       ENDIF
[3182]1848       IF ( bc_radiation_r )  THEN
[151]1849          u_m_r(:,:,:) = u(:,:,nx-1:nx)
1850          v_m_r(:,:,:) = v(:,:,nx-1:nx)
1851          w_m_r(:,:,:) = w(:,:,nx-1:nx)
1852       ENDIF
[3182]1853       IF ( bc_radiation_s )  THEN
[151]1854          u_m_s(:,:,:) = u(:,0:1,:)
1855          v_m_s(:,:,:) = v(:,1:2,:)
1856          w_m_s(:,:,:) = w(:,0:1,:)
1857       ENDIF
[3182]1858       IF ( bc_radiation_n )  THEN
[151]1859          u_m_n(:,:,:) = u(:,ny-1:ny,:)
1860          v_m_n(:,:,:) = v(:,ny-1:ny,:)
1861          w_m_n(:,:,:) = w(:,ny-1:ny,:)
1862       ENDIF
[667]1863       
[151]1864    ENDIF
[680]1865
[667]1866!
1867!-- Calculate the initial volume flow at the right and north boundary
[709]1868    IF ( conserve_volume_flow )  THEN
[151]1869
[767]1870       IF ( use_prescribed_profile_data )  THEN
[667]1871
[1340]1872          volume_flow_initial_l = 0.0_wp
1873          volume_flow_area_l    = 0.0_wp
[732]1874
[667]1875          IF ( nxr == nx )  THEN
1876             DO  j = nys, nyn
[2232]1877                DO  k = nzb+1, nzt
[1788]1878                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
[2232]1879                                              u_init(k) * dzw(k)               &
1880                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1881                                              BTEST( wall_flags_0(k,j,nxr), 1 )&
1882                                            )
1883
1884                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
1885                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1886                                              BTEST( wall_flags_0(k,j,nxr), 1 )&
1887                                            )
[767]1888                ENDDO
1889             ENDDO
1890          ENDIF
1891         
1892          IF ( nyn == ny )  THEN
1893             DO  i = nxl, nxr
[2232]1894                DO  k = nzb+1, nzt
1895                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1896                                              v_init(k) * dzw(k)               &       
1897                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1898                                              BTEST( wall_flags_0(k,nyn,i), 2 )&
1899                                            )
1900                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
1901                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1902                                              BTEST( wall_flags_0(k,nyn,i), 2 )&
1903                                            )
[767]1904                ENDDO
1905             ENDDO
1906          ENDIF
1907
1908#if defined( __parallel )
1909          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1910                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1911          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1912                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1913
1914#else
1915          volume_flow_initial = volume_flow_initial_l
1916          volume_flow_area    = volume_flow_area_l
1917#endif 
1918
1919       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1920
[1340]1921          volume_flow_initial_l = 0.0_wp
1922          volume_flow_area_l    = 0.0_wp
[767]1923
1924          IF ( nxr == nx )  THEN
1925             DO  j = nys, nyn
[2232]1926                DO  k = nzb+1, nzt
[1788]1927                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
[2232]1928                                              hom_sum(k,1,0) * dzw(k)          &
1929                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1930                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
1931                                            )
1932                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
1933                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1934                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
1935                                            )
[667]1936                ENDDO
1937             ENDDO
1938          ENDIF
1939         
1940          IF ( nyn == ny )  THEN
1941             DO  i = nxl, nxr
[2232]1942                DO  k = nzb+1, nzt
[1788]1943                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
[2232]1944                                              hom_sum(k,2,0) * dzw(k)          &       
1945                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1946                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
1947                                            )
1948                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
1949                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1950                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
1951                                            )
[667]1952                ENDDO
1953             ENDDO
1954          ENDIF
1955
[732]1956#if defined( __parallel )
1957          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1958                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1959          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1960                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1961
1962#else
1963          volume_flow_initial = volume_flow_initial_l
1964          volume_flow_area    = volume_flow_area_l
1965#endif 
1966
[667]1967       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1968
[1340]1969          volume_flow_initial_l = 0.0_wp
1970          volume_flow_area_l    = 0.0_wp
[732]1971
[667]1972          IF ( nxr == nx )  THEN
1973             DO  j = nys, nyn
[2232]1974                DO  k = nzb+1, nzt
1975                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1976                                              u(k,j,nx) * dzw(k)               &
1977                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1978                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
1979                                            )
1980                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
1981                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1982                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
1983                                            )
[667]1984                ENDDO
1985             ENDDO
1986          ENDIF
1987         
1988          IF ( nyn == ny )  THEN
1989             DO  i = nxl, nxr
[2232]1990                DO  k = nzb+1, nzt
[1788]1991                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
[2232]1992                                              v(k,ny,i) * dzw(k)               &       
1993                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1994                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
1995                                            )
1996                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
1997                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1998                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
1999                                            )
[667]2000                ENDDO
2001             ENDDO
2002          ENDIF
2003
2004#if defined( __parallel )
[732]2005          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
2006                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
2007          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
2008                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
[667]2009
2010#else
[732]2011          volume_flow_initial = volume_flow_initial_l
2012          volume_flow_area    = volume_flow_area_l
[667]2013#endif 
2014
[732]2015       ENDIF
2016
[151]2017!
[709]2018!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated
2019!--    from u|v_bulk instead
[680]2020       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
2021          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
2022          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
2023       ENDIF
[667]2024
[680]2025    ENDIF
[2232]2026!
[2618]2027!-- Finally, if random_heatflux is set, disturb shf at horizontal
2028!-- surfaces. Actually, this should be done in surface_mod, where all other
2029!-- initializations of surface quantities are done. However, this
2030!-- would create a ring dependency, hence, it is done here. Maybe delete
2031!-- disturb_heatflux and tranfer the respective code directly into the
2032!-- initialization in surface_mod.         
[2232]2033    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
2034         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
[2618]2035 
[2232]2036       IF ( use_surface_fluxes  .AND.  constant_heatflux  .AND.                &
2037            random_heatflux )  THEN
2038          IF ( surf_def_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_def_h(0) )
2039          IF ( surf_lsm_h%ns    >= 1 )  CALL disturb_heatflux( surf_lsm_h    )
2040          IF ( surf_usm_h%ns    >= 1 )  CALL disturb_heatflux( surf_usm_h    )
2041       ENDIF
2042    ENDIF
[680]2043
[787]2044!
[2696]2045!-- Before initializing further modules, compute total sum of active mask
2046!-- grid points and the mean surface level height for each statistic region.
2047!-- ngp_2dh: number of grid points of a horizontal cross section through the
2048!--          total domain
2049!-- ngp_3d:  number of grid points of the total domain
2050    ngp_2dh_outer_l   = 0
2051    ngp_2dh_outer     = 0
2052    ngp_2dh_s_inner_l = 0
2053    ngp_2dh_s_inner   = 0
2054    ngp_2dh_l         = 0
2055    ngp_2dh           = 0
2056    ngp_3d_inner_l    = 0.0_wp
2057    ngp_3d_inner      = 0
2058    ngp_3d            = 0
2059    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
2060
2061    mean_surface_level_height   = 0.0_wp
2062    mean_surface_level_height_l = 0.0_wp
2063!
2064!-- Pre-set masks for regional statistics. Default is the total model domain.
2065!-- Ghost points are excluded because counting values at the ghost boundaries
2066!-- would bias the statistics
2067    rmask = 1.0_wp
2068    rmask(:,nxlg:nxl-1,:) = 0.0_wp;  rmask(:,nxr+1:nxrg,:) = 0.0_wp
2069    rmask(nysg:nys-1,:,:) = 0.0_wp;  rmask(nyn+1:nyng,:,:) = 0.0_wp
[3234]2070
[2696]2071!
[3234]2072!-- Temporary solution to add LSM and radiation time series to the default
2073!-- output
2074    IF ( land_surface  .OR.  radiation )  THEN
2075       IF ( TRIM( radiation_scheme ) == 'rrtmg' )  THEN
2076          dots_num = dots_num + 15
2077       ELSE
2078          dots_num = dots_num + 11
2079       ENDIF
2080    ENDIF
2081!
[2867]2082!-- User-defined initializing actions
2083    CALL user_init
2084!
[2696]2085!-- To do: New concept for these non-topography grid points!
2086    DO  sr = 0, statistic_regions
2087       DO  i = nxl, nxr
2088          DO  j = nys, nyn
2089             IF ( rmask(j,i,sr) == 1.0_wp )  THEN
2090!
2091!--             All xy-grid points
2092                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
2093!
2094!--             Determine mean surface-level height. In case of downward-
2095!--             facing walls are present, more than one surface level exist.
2096!--             In this case, use the lowest surface-level height.
2097                IF ( surf_def_h(0)%start_index(j,i) <=                         &
2098                     surf_def_h(0)%end_index(j,i) )  THEN
2099                   m = surf_def_h(0)%start_index(j,i)
2100                   k = surf_def_h(0)%k(m)
2101                   mean_surface_level_height_l(sr) =                           &
2102                                       mean_surface_level_height_l(sr) + zw(k-1)
2103                ENDIF
2104                IF ( surf_lsm_h%start_index(j,i) <=                            &
2105                     surf_lsm_h%end_index(j,i) )  THEN
2106                   m = surf_lsm_h%start_index(j,i)
2107                   k = surf_lsm_h%k(m)
2108                   mean_surface_level_height_l(sr) =                           &
2109                                       mean_surface_level_height_l(sr) + zw(k-1)
2110                ENDIF
2111                IF ( surf_usm_h%start_index(j,i) <=                            &
2112                     surf_usm_h%end_index(j,i) )  THEN
2113                   m = surf_usm_h%start_index(j,i)
2114                   k = surf_usm_h%k(m)
2115                   mean_surface_level_height_l(sr) =                           &
2116                                       mean_surface_level_height_l(sr) + zw(k-1)
2117                ENDIF
2118
2119                k_surf = k - 1
2120
2121                DO  k = nzb, nzt+1
2122!
2123!--                xy-grid points above topography
2124                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr)     +         &
2125                                  MERGE( 1, 0, BTEST( wall_flags_0(k,j,i), 24 ) )
2126
2127                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) +         &
2128                                  MERGE( 1, 0, BTEST( wall_flags_0(k,j,i), 22 ) )
2129
2130                ENDDO
2131!
2132!--             All grid points of the total domain above topography
2133                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + ( nz - k_surf + 2 )
2134
2135
2136
2137             ENDIF
2138          ENDDO
2139       ENDDO
2140    ENDDO
[2864]2141!
2142!-- Initialize arrays encompassing number of grid-points in inner and outer
2143!-- domains, statistic regions, etc. Mainly used for horizontal averaging
2144!-- of turbulence statistics. Please note, user_init must be called before
2145!-- doing this.   
[2696]2146    sr = statistic_regions + 1
2147#if defined( __parallel )
2148    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2149    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,    &
2150                        comm2d, ierr )
2151    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2152    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr,   &
2153                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
2154    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2155    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),          &
2156                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
2157    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2158    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL,  &
2159                        MPI_SUM, comm2d, ierr )
2160    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
2161    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2162    CALL MPI_ALLREDUCE( mean_surface_level_height_l(0),                        &
2163                        mean_surface_level_height(0), sr, MPI_REAL,            &
2164                        MPI_SUM, comm2d, ierr )
2165    mean_surface_level_height = mean_surface_level_height / REAL( ngp_2dh )
2166#else
2167    ngp_2dh         = ngp_2dh_l
2168    ngp_2dh_outer   = ngp_2dh_outer_l
2169    ngp_2dh_s_inner = ngp_2dh_s_inner_l
2170    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
2171    mean_surface_level_height = mean_surface_level_height_l / REAL( ngp_2dh_l )
2172#endif
2173
2174    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * &
2175             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
2176
2177!
2178!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
2179!-- buoyancy, etc. A zero value will occur for cases where all grid points of
2180!-- the respective subdomain lie below the surface topography
2181    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   ) 
2182    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )),             &
2183                           ngp_3d_inner(:) )
2184    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) 
2185
2186    DEALLOCATE( mean_surface_level_height_l, ngp_2dh_l, ngp_2dh_outer_l,       &
2187                ngp_3d_inner_l, ngp_3d_inner_tmp )
2188
2189!
[3051]2190!-- Initialize nudging if required
2191    IF ( nudging )  CALL nudge_init
2192!
2193!-- Initialize 1D/3D offline-nesting with COSMO model and read data from
2194!-- external file.
[3182]2195    IF ( large_scale_forcing  .OR.  nesting_offline )  CALL lsf_init
[3051]2196!
[2232]2197!-- Initialize surface forcing corresponding to large-scale forcing. Therein,
2198!-- initialize heat-fluxes, etc. via datatype. Revise it later!
2199    IF ( large_scale_forcing .AND. lsf_surf )  THEN
2200       IF ( use_surface_fluxes  .AND.  constant_heatflux )  THEN
2201          CALL ls_forcing_surf ( simulated_time )
2202       ENDIF
2203    ENDIF
2204!
[787]2205!-- Initialize quantities for special advections schemes
2206    CALL init_advec
[680]2207
[667]2208!
[680]2209!-- Impose random perturbation on the horizontal velocity field and then
2210!-- remove the divergences from the velocity field at the initial stage
[1788]2211    IF ( create_disturbances  .AND.  disturbance_energy_limit /= 0.0_wp  .AND. &
2212         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
[680]2213         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
2214
[1402]2215       CALL location_message( 'creating initial disturbances', .FALSE. )
[2232]2216       CALL disturb_field( 'u', tend, u )
2217       CALL disturb_field( 'v', tend, v )
[1402]2218       CALL location_message( 'finished', .TRUE. )
[1384]2219
[1402]2220       CALL location_message( 'calling pressure solver', .FALSE. )
[680]2221       n_sor = nsor_ini
2222       CALL pres
2223       n_sor = nsor
[1402]2224       CALL location_message( 'finished', .TRUE. )
[1384]2225
[680]2226    ENDIF
2227
2228!
[1484]2229!-- If required, initialize quantities needed for the plant canopy model
[2007]2230    IF ( plant_canopy )  THEN
2231       CALL location_message( 'initializing plant canopy model', .FALSE. )   
2232       CALL pcm_init
2233       CALL location_message( 'finished', .TRUE. )
2234    ENDIF
[138]2235
2236!
[1]2237!-- If required, initialize dvrp-software
[1340]2238    IF ( dt_dvrp /= 9999999.9_wp )  CALL init_dvrp
[1]2239
[96]2240    IF ( ocean )  THEN
[1]2241!
[96]2242!--    Initialize quantities needed for the ocean model
2243       CALL init_ocean
[388]2244
[3274]2245    ENDIF
[96]2246!
2247!--    Initialize quantities for handling cloud physics
[849]2248!--    This routine must be called before lpm_init, because
[3274]2249!--    otherwise, array d_exner, needed in data_output_dvrp (called by
[849]2250!--    lpm_init) is not defined.
[3274]2251    IF ( .NOT. ocean )  THEN
2252
2253       ALLOCATE( hyp(nzb:nzt+1) )
2254       ALLOCATE( d_exner(nzb:nzt+1) )
2255       ALLOCATE( exner(nzb:nzt+1) )
2256       ALLOCATE( hyrho(nzb:nzt+1) )
[1849]2257!
[3274]2258!--    Check temperature in case of too large domain height
2259       DO  k = nzb, nzt+1
2260          IF ( ( pt_surface * exner_function(surface_pressure * 100.0_wp) - g/c_p * zu(k) ) < 0.0_wp )  THEN
2261             WRITE( message_string, * )  'absolute temperature < 0.0 at zu(', k, &
2262                                         ') = ', zu(k)
2263             CALL message( 'init_bulk_cloud_model', 'PA0142', 1, 2, 0, 6, 0 )
2264          ENDIF
2265       ENDDO
2266
2267!
2268!--    Calculate vertical profile of the hydrostatic pressure (hyp)
2269       hyp    = barometric_formula(zu, pt_surface * exner_function(surface_pressure * 100.0_wp), surface_pressure * 100.0_wp)
2270       d_exner = exner_function_invers(hyp)
2271       exner = 1.0_wp / exner_function_invers(hyp)
2272       hyrho  = ideal_gas_law_rho_pt(hyp, pt_init)
2273!
2274!--    Compute reference density
2275       rho_surface = ideal_gas_law_rho(surface_pressure * 100.0_wp, pt_surface * exner_function(surface_pressure * 100.0_wp))
2276
[96]2277    ENDIF
[3274]2278!
2279!-- If required, initialize quantities needed for the microphysics module
2280    IF ( bulk_cloud_model )  THEN
2281       CALL bcm_init
2282    ENDIF
[1]2283
2284!
2285!-- If required, initialize particles
[849]2286    IF ( particle_advection )  CALL lpm_init
[1]2287
[1585]2288!
[3159]2289!-- If required, initialize particles
2290    IF ( agents_active )  CALL mas_init
2291
2292!
[1585]2293!-- If required, initialize quantities needed for the LSM
2294    IF ( land_surface )  THEN
2295       CALL location_message( 'initializing land surface model', .FALSE. )
[1817]2296       CALL lsm_init
[1585]2297       CALL location_message( 'finished', .TRUE. )
2298    ENDIF
[1496]2299
[1]2300!
[2696]2301!-- If required, allocate USM and LSM surfaces
2302    IF ( urban_surface )  THEN
2303       CALL location_message( 'initializing and allocating urban surfaces', .FALSE. )
2304       CALL usm_allocate_surface
2305       CALL location_message( 'finished', .TRUE. )
2306    ENDIF
2307!
2308!-- If required, initialize urban surface model
2309    IF ( urban_surface )  THEN
2310       CALL location_message( 'initializing urban surface model', .FALSE. )
2311       CALL usm_init_urban_surface
2312       CALL location_message( 'finished', .TRUE. )
2313    ENDIF
2314
2315!
[1691]2316!-- Initialize surface layer (done after LSM as roughness length are required
2317!-- for initialization
2318    IF ( constant_flux_layer )  THEN
2319       CALL location_message( 'initializing surface layer', .FALSE. )
2320       CALL init_surface_layer_fluxes
2321       CALL location_message( 'finished', .TRUE. )
2322    ENDIF
2323
2324!
[2696]2325!-- If required, set chemical emissions
2326!-- (todo(FK): This should later on be CALLed time-dependently in init_3d_model)
[2977]2327    IF ( air_chemistry )  THEN
[2696]2328       CALL chem_emissions
2329    ENDIF
2330
2331!
2332!-- Initialize radiation processes
[1496]2333    IF ( radiation )  THEN
[2696]2334!
[2977]2335!--    Activate radiation_interactions according to the existence of vertical surfaces and/or trees.
2336!--    The namelist parameter radiation_interactions_on can override this behavior.
2337!--    (This check cannot be performed in check_parameters, because vertical_surfaces_exist is first set in
2338!--    init_surface_arrays.)
2339       IF ( radiation_interactions_on )  THEN
2340          IF ( vertical_surfaces_exist  .OR.  plant_canopy )  THEN
2341             radiation_interactions    = .TRUE.
2342             average_radiation         = .TRUE.
2343          ELSE
2344             radiation_interactions_on = .FALSE.   !< reset namelist parameter: no interactions
2345                                                   !< calculations necessary in case of flat surface
2346          ENDIF
2347       ELSEIF ( vertical_surfaces_exist  .OR.  plant_canopy )  THEN
2348          message_string = 'radiation_interactions_on is set to .FALSE. although '     // &
2349                           'vertical surfaces and/or trees exist. The model will run ' // &
2350                           'without RTM (no shadows, no radiation reflections)'
2351          CALL message( 'init_3d_model', 'PA0348', 0, 1, 0, 6, 0 )
2352       ENDIF
2353!
[2696]2354!--    If required, initialize radiation interactions between surfaces
[2920]2355!--    via sky-view factors. This must be done before radiation is initialized.
[2696]2356       IF ( radiation_interactions )  CALL radiation_interaction_init
2357
2358!
2359!--    Initialize radiation model
[1585]2360       CALL location_message( 'initializing radiation model', .FALSE. )
[1826]2361       CALL radiation_init
[1585]2362       CALL location_message( 'finished', .TRUE. )
[2696]2363
2364!
[2920]2365!--    Find all discretized apparent solar positions for radiation interaction.
2366!--    This must be done after radiation_init.
2367       IF ( radiation_interactions )  CALL radiation_presimulate_solar_pos
2368
2369!
[2696]2370!--    If required, read or calculate and write out the SVF
[2906]2371       IF ( radiation_interactions .AND. read_svf)  THEN
[2696]2372!
2373!--       Read sky-view factors and further required data from file
2374          CALL location_message( '    Start reading SVF from file', .FALSE. )
2375          CALL radiation_read_svf()
2376          CALL location_message( '    Reading SVF from file has finished', .TRUE. )
2377
[2906]2378       ELSEIF ( radiation_interactions .AND. .NOT. read_svf)  THEN
[2696]2379!
2380!--       calculate SFV and CSF
2381          CALL location_message( '    Start calculation of SVF', .FALSE. )
2382          CALL radiation_calc_svf()
2383          CALL location_message( '    Calculation of SVF has finished', .TRUE. )
2384       ENDIF
2385
[2906]2386       IF ( radiation_interactions .AND. write_svf)  THEN
[2696]2387!
2388!--       Write svf, csf svfsurf and csfsurf data to file
[2906]2389          CALL location_message( '    Start writing SVF in file', .FALSE. )
[2696]2390          CALL radiation_write_svf()
[2906]2391          CALL location_message( '    Writing SVF in file has finished', .TRUE. )
[2696]2392       ENDIF
2393
2394!
2395!--    Adjust radiative fluxes. In case of urban and land surfaces, also
2396!--    call an initial interaction.
2397       IF ( radiation_interactions )  THEN
2398          CALL radiation_interaction
2399       ENDIF
[1496]2400    ENDIF
[3234]2401 
[1914]2402!
2403!-- If required, initialize quantities needed for the wind turbine model
2404    IF ( wind_turbine )  THEN
2405       CALL location_message( 'initializing wind turbine model', .FALSE. )
2406       CALL wtm_init
2407       CALL location_message( 'finished', .TRUE. )
2408    ENDIF
[1496]2409
[2817]2410!
2411!-- If required, initialize quantities needed for the gust module
2412    IF ( gust_module_enabled )  THEN
2413       CALL gust_init( dots_label, dots_unit, dots_num, dots_max )
2414    ENDIF
[1914]2415
[1496]2416!
[673]2417!-- Initialize the ws-scheme.   
2418    IF ( ws_scheme_sca .OR. ws_scheme_mom )  CALL ws_init       
[1]2419
2420!
[709]2421!-- Setting weighting factors for calculation of perturbation pressure
[1762]2422!-- and turbulent quantities from the RK substeps
[709]2423    IF ( TRIM(timestep_scheme) == 'runge-kutta-3' )  THEN      ! for RK3-method
2424
[1322]2425       weight_substep(1) = 1._wp/6._wp
2426       weight_substep(2) = 3._wp/10._wp
2427       weight_substep(3) = 8._wp/15._wp
[709]2428
[1322]2429       weight_pres(1)    = 1._wp/3._wp
2430       weight_pres(2)    = 5._wp/12._wp
2431       weight_pres(3)    = 1._wp/4._wp
[709]2432
2433    ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' )  THEN  ! for RK2-method
2434
[1322]2435       weight_substep(1) = 1._wp/2._wp
2436       weight_substep(2) = 1._wp/2._wp
[673]2437         
[1322]2438       weight_pres(1)    = 1._wp/2._wp
2439       weight_pres(2)    = 1._wp/2._wp       
[709]2440
[1001]2441    ELSE                                     ! for Euler-method
[709]2442
[1340]2443       weight_substep(1) = 1.0_wp     
2444       weight_pres(1)    = 1.0_wp                   
[709]2445
[673]2446    ENDIF
2447
2448!
[1]2449!-- Initialize Rayleigh damping factors
[1340]2450    rdf    = 0.0_wp
2451    rdf_sc = 0.0_wp
2452    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
[1788]2453       IF (  .NOT.  ocean )  THEN
[108]2454          DO  k = nzb+1, nzt
2455             IF ( zu(k) >= rayleigh_damping_height )  THEN
[1788]2456                rdf(k) = rayleigh_damping_factor *                             &
[1340]2457                      ( SIN( pi * 0.5_wp * ( zu(k) - rayleigh_damping_height ) &
[1788]2458                             / ( zu(nzt) - rayleigh_damping_height ) )         &
[1]2459                      )**2
[108]2460             ENDIF
2461          ENDDO
2462       ELSE
2463          DO  k = nzt, nzb+1, -1
2464             IF ( zu(k) <= rayleigh_damping_height )  THEN
[1788]2465                rdf(k) = rayleigh_damping_factor *                             &
[1340]2466                      ( SIN( pi * 0.5_wp * ( rayleigh_damping_height - zu(k) ) &
[1788]2467                             / ( rayleigh_damping_height - zu(nzb+1) ) )       &
[108]2468                      )**2
2469             ENDIF
2470          ENDDO
2471       ENDIF
[1]2472    ENDIF
[785]2473    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
[1]2474
2475!
[240]2476!-- Initialize the starting level and the vertical smoothing factor used for
2477!-- the external pressure gradient
[1340]2478    dp_smooth_factor = 1.0_wp
[240]2479    IF ( dp_external )  THEN
2480!
2481!--    Set the starting level dp_level_ind_b only if it has not been set before
2482!--    (e.g. in init_grid).
2483       IF ( dp_level_ind_b == 0 )  THEN
2484          ind_array = MINLOC( ABS( dp_level_b - zu ) )
2485          dp_level_ind_b = ind_array(1) - 1 + nzb 
2486                                        ! MINLOC uses lower array bound 1
2487       ENDIF
2488       IF ( dp_smooth )  THEN
[1340]2489          dp_smooth_factor(:dp_level_ind_b) = 0.0_wp
[240]2490          DO  k = dp_level_ind_b+1, nzt
[1340]2491             dp_smooth_factor(k) = 0.5_wp * ( 1.0_wp + SIN( pi *               &
2492                        ( REAL( k - dp_level_ind_b, KIND=wp ) /                &
2493                          REAL( nzt - dp_level_ind_b, KIND=wp ) - 0.5_wp ) ) )
[240]2494          ENDDO
2495       ENDIF
2496    ENDIF
2497
2498!
[978]2499!-- Initialize damping zone for the potential temperature in case of
2500!-- non-cyclic lateral boundaries. The damping zone has the maximum value
2501!-- at the inflow boundary and decreases to zero at pt_damping_width.
[1340]2502    ptdf_x = 0.0_wp
2503    ptdf_y = 0.0_wp
[1159]2504    IF ( bc_lr_dirrad )  THEN
[996]2505       DO  i = nxl, nxr
[978]2506          IF ( ( i * dx ) < pt_damping_width )  THEN
[1340]2507             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5_wp *              &
2508                            REAL( pt_damping_width - i * dx, KIND=wp ) / (     &
[1788]2509                            REAL( pt_damping_width, KIND=wp ) ) ) )**2 
[73]2510          ENDIF
2511       ENDDO
[1159]2512    ELSEIF ( bc_lr_raddir )  THEN
[996]2513       DO  i = nxl, nxr
[978]2514          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
[1322]2515             ptdf_x(i) = pt_damping_factor *                                   &
[1340]2516                         SIN( pi * 0.5_wp *                                    &
2517                                 ( ( i - nx ) * dx + pt_damping_width ) /      &
2518                                 REAL( pt_damping_width, KIND=wp ) )**2
[73]2519          ENDIF
[978]2520       ENDDO 
[1159]2521    ELSEIF ( bc_ns_dirrad )  THEN
[996]2522       DO  j = nys, nyn
[978]2523          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
[1322]2524             ptdf_y(j) = pt_damping_factor *                                   &
[1340]2525                         SIN( pi * 0.5_wp *                                    &
2526                                 ( ( j - ny ) * dy + pt_damping_width ) /      &
2527                                 REAL( pt_damping_width, KIND=wp ) )**2
[1]2528          ENDIF
[978]2529       ENDDO 
[1159]2530    ELSEIF ( bc_ns_raddir )  THEN
[996]2531       DO  j = nys, nyn
[978]2532          IF ( ( j * dy ) < pt_damping_width )  THEN
[1322]2533             ptdf_y(j) = pt_damping_factor *                                   &
[1340]2534                         SIN( pi * 0.5_wp *                                    &
2535                                ( pt_damping_width - j * dy ) /                &
2536                                REAL( pt_damping_width, KIND=wp ) )**2
[1]2537          ENDIF
[73]2538       ENDDO
[1]2539    ENDIF
2540!
[2864]2541!-- Check if maximum number of allowed timeseries is exceeded
[51]2542    IF ( dots_num > dots_max )  THEN
[1788]2543       WRITE( message_string, * ) 'number of time series quantities exceeds',  &
2544                                  ' its maximum of dots_max = ', dots_max,     &
[3046]2545                                  '&Please increase dots_max in modules.f90.'
[254]2546       CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 )   
[51]2547    ENDIF
2548
[1]2549!
2550!-- Input binary data file is not needed anymore. This line must be placed
2551!-- after call of user_init!
2552    CALL close_file( 13 )
[2934]2553!
2554!-- In case of nesting, put an barrier to assure that all parent and child
2555!-- domains finished initialization.
2556#if defined( __parallel )
2557    IF ( nested_run )  CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
2558#endif
[1]2559
[2934]2560
[1402]2561    CALL location_message( 'leaving init_3d_model', .TRUE. )
[1]2562
2563 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.