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

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

Move initialization of STG behind initialization of offline nesting; Bugfix in STG in case of very early restart; calculation of scaling parameters used for parametrization of synthetic turbulence profiles improved; in offline nesting, set boundary value at upper-left and upper-south grid points for u- and v-component, respectively

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