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

Last change on this file since 3877 was 3849, checked in by knoop, 5 years ago

Bugfix: added proper OpenACC support to disturb_field

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