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

Last change on this file since 4039 was 4028, checked in by schwenkel, 6 years ago

Further modularization of particle code components

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