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

Last change on this file since 3458 was 3458, checked in by kanani, 5 years ago

Reintegrated fixes/changes from branch chemistry

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