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

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

restructure/add location/debug messages

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