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

Last change on this file since 4130 was 4130, checked in by suehring, 5 years ago

Effectively reduce 3D initialization from dynamic driver to 1D initial profiles. This is because 3D initialization produces structures in the w-component that are correlated with the processor grid for some unknown reason

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