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

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

Consider time-dependent geostrophic wind components in offline nesting

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