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

Last change on this file since 4055 was 4048, checked in by knoop, 5 years ago

Moved turbulence_closure_mod calls into module_interface

  • 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.0 KB
Line 
1!> @file init_3d_model.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2019 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: init_3d_model.f90 4048 2019-06-21 21:00:21Z 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    IMPLICIT NONE
697
698    INTEGER(iwp) ::  i             !< grid index in x direction
699    INTEGER(iwp) ::  ind_array(1)  !< dummy used to determine start index for external pressure forcing
700    INTEGER(iwp) ::  j             !< grid index in y direction
701    INTEGER(iwp) ::  k             !< grid index in z direction
702    INTEGER(iwp) ::  k_surf        !< surface level index
703    INTEGER(iwp) ::  m             !< index of surface element in surface data type
704    INTEGER(iwp) ::  sr            !< index of statistic region
705
706    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  ngp_2dh_l  !< toal number of horizontal grid points in statistical region on subdomain
707
708    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l    !< number of horizontal non-wall bounded grid points on subdomain
709    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner_l  !< number of horizontal non-topography grid points on subdomain
710
711    REAL(wp)     ::  t_surface !< air temperature at the surface
712
713    REAL(wp), DIMENSION(:), ALLOCATABLE ::  init_l        !< dummy array used for averaging 3D data to obtain inital profiles
714    REAL(wp), DIMENSION(:), ALLOCATABLE ::  p_hydrostatic !< hydrostatic pressure
715
716    INTEGER(iwp) ::  l       !< loop variable
717    INTEGER(iwp) ::  nzt_l   !< index of top PE boundary for multigrid level
718    REAL(wp) ::  dx_l !< grid spacing along x on different multigrid level
719    REAL(wp) ::  dy_l !< grid spacing along y on different multigrid level
720
721    REAL(wp), DIMENSION(1:3) ::  volume_flow_area_l     !< area of lateral and top model domain surface on local subdomain
722    REAL(wp), DIMENSION(1:3) ::  volume_flow_initial_l  !< initial volume flow into model domain
723
724    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height_l !< mean surface level height on subdomain
725    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l    !< total number of non-topography grid points on subdomain
726    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_tmp  !< total number of non-topography grid points
727
728    INTEGER(iwp) ::  nz_u_shift   !< topography-top index on u-grid, used to vertically shift initial profiles
729    INTEGER(iwp) ::  nz_v_shift   !< topography-top index on v-grid, used to vertically shift initial profiles
730    INTEGER(iwp) ::  nz_w_shift   !< topography-top index on w-grid, used to vertically shift initial profiles
731    INTEGER(iwp) ::  nz_s_shift   !< topography-top index on scalar-grid, used to vertically shift initial profiles
732    INTEGER(iwp) ::  nz_u_shift_l !< topography-top index on u-grid, used to vertically shift initial profiles
733    INTEGER(iwp) ::  nz_v_shift_l !< topography-top index on v-grid, used to vertically shift initial profiles
734    INTEGER(iwp) ::  nz_w_shift_l !< topography-top index on w-grid, used to vertically shift initial profiles
735    INTEGER(iwp) ::  nz_s_shift_l !< topography-top index on scalar-grid, used to vertically shift initial profiles
736
737
738    CALL location_message( 'model initialization', 'start' )
739
740    IF ( debug_output )  CALL debug_message( 'allocating arrays', 'start' )
741!
742!-- Allocate arrays
743    ALLOCATE( mean_surface_level_height(0:statistic_regions),                  &
744              mean_surface_level_height_l(0:statistic_regions),                &
745              ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions),    &
746              ngp_3d(0:statistic_regions),                                     &
747              ngp_3d_inner(0:statistic_regions),                               &
748              ngp_3d_inner_l(0:statistic_regions),                             &
749              ngp_3d_inner_tmp(0:statistic_regions),                           &
750              sums_divnew_l(0:statistic_regions),                              &
751              sums_divold_l(0:statistic_regions) )
752    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
753    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                    &
754              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),                  &
755              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),                  &
756              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),                &
757              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),                  &
758              sums(nzb:nzt+1,pr_palm+max_pr_user+max_pr_cs),                   &
759              sums_l(nzb:nzt+1,pr_palm+max_pr_user+max_pr_cs,0:threads_per_task-1),      &
760              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1),    &
761              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions) )
762    ALLOCATE( ts_value(dots_max,0:statistic_regions) )
763    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
764
765    ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr),                                    &
766              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
767              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
768
769    ALLOCATE( pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
770              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
771              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
772              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
773              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
774              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
775              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
776              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
777              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
778              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
779              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
780    IF (  .NOT.  neutral )  THEN
781       ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
782    ENDIF
783!
784!-- Pre-set masks for regional statistics. Default is the total model domain.
785!-- Ghost points are excluded because counting values at the ghost boundaries
786!-- would bias the statistics
787    rmask = 1.0_wp
788    rmask(:,nxlg:nxl-1,:) = 0.0_wp;  rmask(:,nxr+1:nxrg,:) = 0.0_wp
789    rmask(nysg:nys-1,:,:) = 0.0_wp;  rmask(nyn+1:nyng,:,:) = 0.0_wp
790!
791!-- Following array is required for perturbation pressure within the iterative
792!-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds
793!-- the weighted average of the substeps and cannot be used in the Poisson
794!-- solver.
795    IF ( psolver == 'sor' )  THEN
796       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
797    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
798!
799!--    For performance reasons, multigrid is using one ghost layer only
800       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
801    ENDIF
802
803!
804!-- Array for storing constant coeffficients of the tridiagonal solver
805    IF ( psolver == 'poisfft' )  THEN
806       ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) )
807       ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) )
808    ENDIF
809
810    IF ( humidity )  THEN
811!
812!--    3D-humidity
813       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
814                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
815                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
816                 vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) 
817    ENDIF   
818   
819    IF ( passive_scalar )  THEN
820
821!
822!--    3D scalar arrays
823       ALLOCATE( s_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
824                 s_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
825                 s_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
826
827    ENDIF
828
829!
830!-- Allocate and set 1d-profiles for Stokes drift velocity. It may be set to
831!-- non-zero values later in ocean_init
832    ALLOCATE( u_stokes_zu(nzb:nzt+1), u_stokes_zw(nzb:nzt+1),                  &
833              v_stokes_zu(nzb:nzt+1), v_stokes_zw(nzb:nzt+1) )
834    u_stokes_zu(:) = 0.0_wp
835    u_stokes_zw(:) = 0.0_wp
836    v_stokes_zu(:) = 0.0_wp
837    v_stokes_zw(:) = 0.0_wp
838
839!
840!-- Allocation of anelastic and Boussinesq approximation specific arrays
841    ALLOCATE( p_hydrostatic(nzb:nzt+1) )
842    ALLOCATE( rho_air(nzb:nzt+1) )
843    ALLOCATE( rho_air_zw(nzb:nzt+1) )
844    ALLOCATE( drho_air(nzb:nzt+1) )
845    ALLOCATE( drho_air_zw(nzb:nzt+1) )
846!
847!-- Density profile calculation for anelastic approximation
848    t_surface = pt_surface * ( surface_pressure / 1000.0_wp )**( r_d / c_p )
849    IF ( TRIM( approximation ) == 'anelastic' ) THEN
850       DO  k = nzb, nzt+1
851          p_hydrostatic(k)    = surface_pressure * 100.0_wp *                  &
852                                ( 1 - ( g * zu(k) ) / ( c_p * t_surface )      &
853                                )**( c_p / r_d )
854          rho_air(k)          = ( p_hydrostatic(k) *                           &
855                                  ( 100000.0_wp / p_hydrostatic(k)             &
856                                  )**( r_d / c_p )                             &
857                                ) / ( r_d * pt_init(k) )
858       ENDDO
859       DO  k = nzb, nzt
860          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
861       ENDDO
862       rho_air_zw(nzt+1)  = rho_air_zw(nzt)                                    &
863                            + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
864    ELSE
865       DO  k = nzb, nzt+1
866          p_hydrostatic(k)    = surface_pressure * 100.0_wp *                  &
867                                ( 1 - ( g * zu(nzb) ) / ( c_p * t_surface )    &
868                                )**( c_p / r_d )
869          rho_air(k)          = ( p_hydrostatic(k) *                           &
870                                  ( 100000.0_wp / p_hydrostatic(k)             &
871                                  )**( r_d / c_p )                             &
872                                ) / ( r_d * pt_init(nzb) )
873       ENDDO
874       DO  k = nzb, nzt
875          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
876       ENDDO
877       rho_air_zw(nzt+1)  = rho_air_zw(nzt)                                    &
878                            + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
879    ENDIF
880!
881!-- compute the inverse density array in order to avoid expencive divisions
882    drho_air    = 1.0_wp / rho_air
883    drho_air_zw = 1.0_wp / rho_air_zw
884
885!
886!-- Allocation of flux conversion arrays
887    ALLOCATE( heatflux_input_conversion(nzb:nzt+1) )
888    ALLOCATE( waterflux_input_conversion(nzb:nzt+1) )
889    ALLOCATE( momentumflux_input_conversion(nzb:nzt+1) )
890    ALLOCATE( heatflux_output_conversion(nzb:nzt+1) )
891    ALLOCATE( waterflux_output_conversion(nzb:nzt+1) )
892    ALLOCATE( momentumflux_output_conversion(nzb:nzt+1) )
893
894!
895!-- calculate flux conversion factors according to approximation and in-/output mode
896    DO  k = nzb, nzt+1
897
898        IF ( TRIM( flux_input_mode ) == 'kinematic' )  THEN
899            heatflux_input_conversion(k)      = rho_air_zw(k)
900            waterflux_input_conversion(k)     = rho_air_zw(k)
901            momentumflux_input_conversion(k)  = rho_air_zw(k)
902        ELSEIF ( TRIM( flux_input_mode ) == 'dynamic' ) THEN
903            heatflux_input_conversion(k)      = 1.0_wp / c_p
904            waterflux_input_conversion(k)     = 1.0_wp / l_v
905            momentumflux_input_conversion(k)  = 1.0_wp
906        ENDIF
907
908        IF ( TRIM( flux_output_mode ) == 'kinematic' )  THEN
909            heatflux_output_conversion(k)     = drho_air_zw(k)
910            waterflux_output_conversion(k)    = drho_air_zw(k)
911            momentumflux_output_conversion(k) = drho_air_zw(k)
912        ELSEIF ( TRIM( flux_output_mode ) == 'dynamic' ) THEN
913            heatflux_output_conversion(k)     = c_p
914            waterflux_output_conversion(k)    = l_v
915            momentumflux_output_conversion(k) = 1.0_wp
916        ENDIF
917
918        IF ( .NOT. humidity ) THEN
919            waterflux_input_conversion(k)  = 1.0_wp
920            waterflux_output_conversion(k) = 1.0_wp
921        ENDIF
922
923    ENDDO
924
925!
926!-- In case of multigrid method, compute grid lengths and grid factors for the
927!-- grid levels with respective density on each grid
928    IF ( psolver(1:9) == 'multigrid' )  THEN
929
930       ALLOCATE( ddx2_mg(maximum_grid_level) )
931       ALLOCATE( ddy2_mg(maximum_grid_level) )
932       ALLOCATE( dzu_mg(nzb+1:nzt+1,maximum_grid_level) )
933       ALLOCATE( dzw_mg(nzb+1:nzt+1,maximum_grid_level) )
934       ALLOCATE( f1_mg(nzb+1:nzt,maximum_grid_level) )
935       ALLOCATE( f2_mg(nzb+1:nzt,maximum_grid_level) )
936       ALLOCATE( f3_mg(nzb+1:nzt,maximum_grid_level) )
937       ALLOCATE( rho_air_mg(nzb:nzt+1,maximum_grid_level) )
938       ALLOCATE( rho_air_zw_mg(nzb:nzt+1,maximum_grid_level) )
939
940       dzu_mg(:,maximum_grid_level) = dzu
941       rho_air_mg(:,maximum_grid_level) = rho_air
942!       
943!--    Next line to ensure an equally spaced grid.
944       dzu_mg(1,maximum_grid_level) = dzu(2)
945       rho_air_mg(nzb,maximum_grid_level) = rho_air(nzb) +                     &
946                                             (rho_air(nzb) - rho_air(nzb+1))
947
948       dzw_mg(:,maximum_grid_level) = dzw
949       rho_air_zw_mg(:,maximum_grid_level) = rho_air_zw
950       nzt_l = nzt
951       DO  l = maximum_grid_level-1, 1, -1
952           dzu_mg(nzb+1,l) = 2.0_wp * dzu_mg(nzb+1,l+1)
953           dzw_mg(nzb+1,l) = 2.0_wp * dzw_mg(nzb+1,l+1)
954           rho_air_mg(nzb,l)    = rho_air_mg(nzb,l+1) + (rho_air_mg(nzb,l+1) - rho_air_mg(nzb+1,l+1))
955           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))
956           rho_air_mg(nzb+1,l)    = rho_air_mg(nzb+1,l+1)
957           rho_air_zw_mg(nzb+1,l) = rho_air_zw_mg(nzb+1,l+1)
958           nzt_l = nzt_l / 2
959           DO  k = 2, nzt_l+1
960              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
961              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
962              rho_air_mg(k,l)    = rho_air_mg(2*k-1,l+1)
963              rho_air_zw_mg(k,l) = rho_air_zw_mg(2*k-1,l+1)
964           ENDDO
965       ENDDO
966
967       nzt_l = nzt
968       dx_l  = dx
969       dy_l  = dy
970       DO  l = maximum_grid_level, 1, -1
971          ddx2_mg(l) = 1.0_wp / dx_l**2
972          ddy2_mg(l) = 1.0_wp / dy_l**2
973          DO  k = nzb+1, nzt_l
974             f2_mg(k,l) = rho_air_zw_mg(k,l) / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
975             f3_mg(k,l) = rho_air_zw_mg(k-1,l) / ( dzu_mg(k,l)   * dzw_mg(k,l) )
976             f1_mg(k,l) = 2.0_wp * ( ddx2_mg(l) + ddy2_mg(l) ) &
977                          * rho_air_mg(k,l) + f2_mg(k,l) + f3_mg(k,l)
978          ENDDO
979          nzt_l = nzt_l / 2
980          dx_l  = dx_l * 2.0_wp
981          dy_l  = dy_l * 2.0_wp
982       ENDDO
983
984    ENDIF
985
986!
987!-- 1D-array for large scale subsidence velocity
988    IF ( .NOT. ALLOCATED( w_subs ) )  THEN
989       ALLOCATE ( w_subs(nzb:nzt+1) )
990       w_subs = 0.0_wp
991    ENDIF
992
993!
994!-- Arrays to store velocity data from t-dt and the phase speeds which
995!-- are needed for radiation boundary conditions
996    IF ( bc_radiation_l )  THEN
997       ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2),                               &
998                 v_m_l(nzb:nzt+1,nysg:nyng,0:1),                               &
999                 w_m_l(nzb:nzt+1,nysg:nyng,0:1) )
1000    ENDIF
1001    IF ( bc_radiation_r )  THEN
1002       ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
1003                 v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
1004                 w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) )
1005    ENDIF
1006    IF ( bc_radiation_l  .OR.  bc_radiation_r )  THEN
1007       ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng),           &
1008                 c_w(nzb:nzt+1,nysg:nyng) )
1009    ENDIF
1010    IF ( bc_radiation_s )  THEN
1011       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg),                               &
1012                 v_m_s(nzb:nzt+1,1:2,nxlg:nxrg),                               &
1013                 w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) )
1014    ENDIF
1015    IF ( bc_radiation_n )  THEN
1016       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
1017                 v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
1018                 w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) )
1019    ENDIF
1020    IF ( bc_radiation_s  .OR.  bc_radiation_n )  THEN
1021       ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg),           &
1022                 c_w(nzb:nzt+1,nxlg:nxrg) )
1023    ENDIF
1024    IF ( bc_radiation_l  .OR.  bc_radiation_r  .OR.  bc_radiation_s  .OR.      &
1025         bc_radiation_n )  THEN
1026       ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) )                   
1027       ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) )
1028    ENDIF
1029
1030!
1031!-- Initial assignment of the pointers
1032    IF ( .NOT. neutral )  THEN
1033       pt => pt_1;  pt_p => pt_2;  tpt_m => pt_3
1034    ELSE
1035       pt => pt_1;  pt_p => pt_1;  tpt_m => pt_3
1036    ENDIF
1037    u  => u_1;   u_p  => u_2;   tu_m  => u_3
1038    v  => v_1;   v_p  => v_2;   tv_m  => v_3
1039    w  => w_1;   w_p  => w_2;   tw_m  => w_3
1040
1041    IF ( humidity )  THEN
1042       q => q_1;  q_p => q_2;  tq_m => q_3
1043       vpt  => vpt_1
1044    ENDIF
1045   
1046    IF ( passive_scalar )  THEN
1047       s => s_1;  s_p => s_2;  ts_m => s_3
1048    ENDIF   
1049
1050!
1051!-- Initialize surface arrays
1052    CALL init_surface_arrays
1053!
1054!-- Allocate arrays for other modules
1055    CALL module_interface_init_arrays
1056
1057
1058!
1059!-- Allocate arrays containing the RK coefficient for calculation of
1060!-- perturbation pressure and turbulent fluxes. At this point values are
1061!-- set for pressure calculation during initialization (where no timestep
1062!-- is done). Further below the values needed within the timestep scheme
1063!-- will be set.
1064    ALLOCATE( weight_substep(1:intermediate_timestep_count_max),               &
1065              weight_pres(1:intermediate_timestep_count_max) )
1066    weight_substep = 1.0_wp
1067    weight_pres    = 1.0_wp
1068    intermediate_timestep_count = 0  ! needed when simulated_time = 0.0
1069       
1070    IF ( debug_output )  CALL debug_message( 'allocating arrays', 'end' )
1071
1072!
1073!-- Initialize time series
1074    ts_value = 0.0_wp
1075
1076!
1077!-- Initialize local summation arrays for routine flow_statistics.
1078!-- This is necessary because they may not yet have been initialized when they
1079!-- are called from flow_statistics (or - depending on the chosen model run -
1080!-- are never initialized)
1081    sums_divnew_l      = 0.0_wp
1082    sums_divold_l      = 0.0_wp
1083    sums_l_l           = 0.0_wp
1084    sums_wsts_bc_l     = 0.0_wp
1085   
1086!
1087!-- Initialize model variables
1088    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
1089         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
1090!
1091!--    Initialization with provided input data derived from larger-scale model
1092       IF ( INDEX( initializing_actions, 'inifor' ) /= 0 )  THEN
1093          IF ( debug_output )  CALL debug_message( 'initializing with INIFOR', 'start' )
1094!
1095!--       Read initial 1D profiles or 3D data from NetCDF file, depending
1096!--       on the provided level-of-detail.
1097!--       At the moment, only u, v, w, pt and q are provided.
1098          CALL netcdf_data_input_init_3d
1099!
1100!--       Please note, Inifor provides data from nzb+1 to nzt.
1101!--       Bottom and top boundary conditions for Inifor profiles are already
1102!--       set (just after reading), so that this is not necessary here.
1103!--       Depending on the provided level-of-detail, initial Inifor data is
1104!--       either stored on data type (lod=1), or directly on 3D arrays (lod=2).
1105!--       In order to obtain also initial profiles in case of lod=2 (which
1106!--       is required for e.g. damping), average over 3D data.
1107          IF( init_3d%lod_u == 1 )  THEN
1108             u_init = init_3d%u_init
1109          ELSEIF( init_3d%lod_u == 2 )  THEN
1110             ALLOCATE( init_l(nzb:nzt+1) ) 
1111             DO  k = nzb, nzt+1
1112                init_l(k) = SUM( u(k,nys:nyn,nxl:nxr) )
1113             ENDDO
1114             init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
1115
1116#if defined( __parallel )
1117             CALL MPI_ALLREDUCE( init_l, u_init, nzt+1-nzb+1,                  &
1118                                 MPI_REAL, MPI_SUM, comm2d, ierr )
1119#else
1120             u_init = init_l
1121#endif
1122             DEALLOCATE( init_l )
1123
1124          ENDIF
1125           
1126          IF( init_3d%lod_v == 1 )  THEN 
1127             v_init = init_3d%v_init
1128          ELSEIF( init_3d%lod_v == 2 )  THEN
1129             ALLOCATE( init_l(nzb:nzt+1) ) 
1130             DO  k = nzb, nzt+1
1131                init_l(k) = SUM( v(k,nys:nyn,nxl:nxr) )
1132             ENDDO
1133             init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
1134
1135#if defined( __parallel )
1136             CALL MPI_ALLREDUCE( init_l, v_init, nzt+1-nzb+1,                  &
1137                                 MPI_REAL, MPI_SUM, comm2d, ierr )
1138#else
1139             v_init = init_l
1140#endif
1141             DEALLOCATE( init_l )
1142          ENDIF
1143          IF( .NOT. neutral )  THEN
1144             IF( init_3d%lod_pt == 1 )  THEN
1145                pt_init = init_3d%pt_init
1146             ELSEIF( init_3d%lod_pt == 2 )  THEN
1147                ALLOCATE( init_l(nzb:nzt+1) ) 
1148                DO  k = nzb, nzt+1
1149                   init_l(k) = SUM( pt(k,nys:nyn,nxl:nxr) )
1150                ENDDO
1151                init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
1152
1153#if defined( __parallel )
1154                CALL MPI_ALLREDUCE( init_l, pt_init, nzt+1-nzb+1,               &
1155                                    MPI_REAL, MPI_SUM, comm2d, ierr )
1156#else
1157                pt_init = init_l
1158#endif
1159                DEALLOCATE( init_l )
1160             ENDIF
1161          ENDIF
1162
1163
1164          IF( humidity )  THEN
1165             IF( init_3d%lod_q == 1 )  THEN
1166                q_init = init_3d%q_init
1167             ELSEIF( init_3d%lod_q == 2 )  THEN
1168                ALLOCATE( init_l(nzb:nzt+1) ) 
1169                DO  k = nzb, nzt+1
1170                   init_l(k) = SUM( q(k,nys:nyn,nxl:nxr) )
1171                ENDDO
1172                init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
1173
1174#if defined( __parallel )
1175                CALL MPI_ALLREDUCE( init_l, q_init, nzt+1-nzb+1,               &
1176                                    MPI_REAL, MPI_SUM, comm2d, ierr )
1177#else
1178                q_init = init_l
1179#endif
1180                DEALLOCATE( init_l )
1181             ENDIF
1182          ENDIF
1183
1184!
1185!--       Write initial profiles onto 3D arrays. Note, only in case of lod = 1,
1186!--       for lod = 2 data is already on 3D arrays.   
1187          DO  i = nxlg, nxrg
1188             DO  j = nysg, nyng
1189                IF( init_3d%lod_u == 1 )  u(:,j,i) = u_init(:)
1190                IF( init_3d%lod_v == 1 )  v(:,j,i) = v_init(:)
1191                IF( .NOT. neutral  .AND.  init_3d%lod_pt == 1 )                &
1192                   pt(:,j,i) = pt_init(:)
1193                IF( humidity  .AND.  init_3d%lod_q == 1 )  q(:,j,i) = q_init(:)
1194             ENDDO
1195          ENDDO
1196!
1197!--       Exchange ghost points and set boundary conditions in case of
1198!--       level-of-detail = 2
1199          IF( init_3d%lod_u == 2 )  CALL exchange_horiz( u, nbgp )
1200          IF( init_3d%lod_v == 2 )  CALL exchange_horiz( v, nbgp )
1201          IF( init_3d%lod_w == 2 )  CALL exchange_horiz( w, nbgp )
1202          IF( .NOT. neutral  .AND.  init_3d%lod_pt == 2 )                      &
1203             CALL exchange_horiz( pt, nbgp )
1204          IF( humidity  .AND.  init_3d%lod_q == 2 )                            &
1205             CALL exchange_horiz( q, nbgp )
1206         
1207          IF ( bc_dirichlet_l )  THEN
1208             DO  j = nysg, nyng
1209                DO  k = nzb, nzt+1
1210                   IF( init_3d%lod_u == 2 )  u(k,j,nxlg:nxl)   = u(k,j,nxlu)
1211                   IF( init_3d%lod_v == 2 )  v(k,j,nxlg:nxl-1) = v(k,j,nxl)
1212                   IF( init_3d%lod_w == 2 )  w(k,j,nxlg:nxl-1) = w(k,j,nxl)
1213                   IF( .NOT. neutral  .AND.  init_3d%lod_pt == 2 )                      &
1214                      pt(k,j,nxlg:nxl-1) = pt(k,j,nxl)
1215                   IF( humidity  .AND.  init_3d%lod_q == 2 )                            &
1216                      q(k,j,nxlg:nxl-1)  = q(k,j,nxl) 
1217                ENDDO
1218             ENDDO
1219          ENDIF
1220          IF ( bc_dirichlet_r )  THEN
1221             DO  j = nysg, nyng
1222                DO  k = nzb, nzt+1
1223                   IF( init_3d%lod_u == 2 )  u(k,j,nxr+1:nxrg) = u(k,j,nxr)
1224                   IF( init_3d%lod_v == 2 )  v(k,j,nxr+1:nxrg) = v(k,j,nxr)
1225                   IF( init_3d%lod_w == 2 )  w(k,j,nxr+1:nxrg) = w(k,j,nxr)
1226                   IF( .NOT. neutral  .AND.  init_3d%lod_pt == 2 )                      &
1227                      pt(k,j,nxr+1:nxrg) = pt(k,j,nxr)
1228                   IF( humidity  .AND.  init_3d%lod_q == 2 )                            &
1229                      q(k,j,nxr+1:nxrg)  = q(k,j,nxr) 
1230                ENDDO
1231             ENDDO
1232          ENDIF
1233          IF ( bc_dirichlet_s )  THEN
1234             DO  i = nxlg, nxrg
1235                DO  k = nzb, nzt+1
1236                   IF( init_3d%lod_u == 2 )  u(k,nysg:nys-1,i) = u(k,nys,i)
1237                   IF( init_3d%lod_v == 2 )  v(k,nysg:nys,i)   = v(k,nysv,i)
1238                   IF( init_3d%lod_w == 2 )  w(k,nysg:nys-1,i) = w(k,nys,i)
1239                   IF( .NOT. neutral  .AND.  init_3d%lod_pt == 2 )                      &
1240                      pt(k,nysg:nys-1,i) = pt(k,nys,i)
1241                   IF( humidity  .AND.  init_3d%lod_q == 2 )                            &
1242                      q(k,nysg:nys-1,i)  = q(k,nys,i)
1243                ENDDO
1244             ENDDO
1245          ENDIF
1246          IF ( bc_dirichlet_n )  THEN
1247             DO  i = nxlg, nxrg
1248                DO  k = nzb, nzt+1
1249                   IF( init_3d%lod_u == 2 )  u(k,nyn+1:nyng,i) = u(k,nyn,i)
1250                   IF( init_3d%lod_v == 2 )  v(k,nyn+1:nyng,i) = v(k,nyn,i)
1251                   IF( init_3d%lod_w == 2 )  w(k,nyn+1:nyng,i) = w(k,nyn,i)
1252                   IF( .NOT. neutral  .AND.  init_3d%lod_pt == 2 )                      &
1253                      pt(k,nyn+1:nyng,i) = pt(k,nyn,i)
1254                   IF( humidity  .AND.  init_3d%lod_q == 2 )                            &
1255                      q(k,nyn+1:nyng,i)  = q(k,nyn,i)
1256                ENDDO
1257             ENDDO
1258          ENDIF
1259!
1260!--       Set geostrophic wind components. 
1261          IF ( init_3d%from_file_ug )  THEN
1262             ug(:) = init_3d%ug_init(:)
1263          ENDIF
1264          IF ( init_3d%from_file_vg )  THEN
1265             vg(:) = init_3d%vg_init(:)
1266          ENDIF
1267!
1268!--       Set bottom and top boundary condition for geostrophic wind
1269          ug(nzt+1) = ug(nzt)
1270          vg(nzt+1) = vg(nzt)
1271          ug(nzb)   = ug(nzb+1)
1272          vg(nzb)   = vg(nzb+1)
1273!
1274!--       Set inital w to 0
1275          w = 0.0_wp
1276
1277          IF ( passive_scalar )  THEN
1278             DO  i = nxlg, nxrg
1279                DO  j = nysg, nyng
1280                   s(:,j,i) = s_init
1281                ENDDO
1282             ENDDO
1283          ENDIF
1284
1285!
1286!--       Set velocity components at non-atmospheric / oceanic grid points to
1287!--       zero.
1288          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1289          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )
1290          w = MERGE( w, 0.0_wp, BTEST( wall_flags_0, 3 ) )
1291!
1292!--       Initialize surface variables, e.g. friction velocity, momentum
1293!--       fluxes, etc.
1294          CALL init_surfaces
1295
1296          IF ( debug_output )  CALL debug_message( 'initializing with INIFOR', 'end' )
1297!
1298!--    Initialization via computed 1D-model profiles
1299       ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1300
1301          IF ( debug_output )  CALL debug_message( 'initializing with 1D model profiles', 'start' )
1302!
1303!--       Use solutions of the 1D model as initial profiles,
1304!--       start 1D model
1305          CALL init_1d_model
1306!
1307!--       Transfer initial profiles to the arrays of the 3D model
1308          DO  i = nxlg, nxrg
1309             DO  j = nysg, nyng
1310                pt(:,j,i) = pt_init
1311                u(:,j,i)  = u1d
1312                v(:,j,i)  = v1d
1313             ENDDO
1314          ENDDO
1315
1316          IF ( humidity )  THEN
1317             DO  i = nxlg, nxrg
1318                DO  j = nysg, nyng
1319                   q(:,j,i) = q_init
1320                ENDDO
1321             ENDDO
1322          ENDIF
1323
1324          IF ( passive_scalar )  THEN
1325             DO  i = nxlg, nxrg
1326                DO  j = nysg, nyng
1327                   s(:,j,i) = s_init
1328                ENDDO
1329             ENDDO   
1330          ENDIF
1331!
1332!--          Store initial profiles for output purposes etc.
1333          IF ( .NOT. constant_diffusion )  THEN
1334             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
1335          ENDIF
1336!
1337!--       Set velocities back to zero
1338          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1339          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )         
1340!
1341!--       WARNING: The extra boundary conditions set after running the
1342!--       -------  1D model impose an error on the divergence one layer
1343!--                below the topography; need to correct later
1344!--       ATTENTION: Provisional correction for Piacsek & Williams
1345!--       ---------  advection scheme: keep u and v zero one layer below
1346!--                  the topography.
1347          IF ( ibc_uv_b == 1 )  THEN
1348!
1349!--          Neumann condition
1350             DO  i = nxl-1, nxr+1
1351                DO  j = nys-1, nyn+1
1352                   u(nzb,j,i) = u(nzb+1,j,i)
1353                   v(nzb,j,i) = v(nzb+1,j,i)
1354                ENDDO
1355             ENDDO
1356
1357          ENDIF
1358!
1359!--       Initialize surface variables, e.g. friction velocity, momentum
1360!--       fluxes, etc.
1361          CALL init_surfaces
1362
1363          IF ( debug_output )  CALL debug_message( 'initializing with 1D model profiles', 'end' )
1364
1365       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )    &
1366       THEN
1367
1368          IF ( debug_output )  CALL debug_message( 'initializing with constant profiles', 'start' )
1369
1370!
1371!--       Use constructed initial profiles (velocity constant with height,
1372!--       temperature profile with constant gradient)
1373          DO  i = nxlg, nxrg
1374             DO  j = nysg, nyng
1375                pt(:,j,i) = pt_init
1376                u(:,j,i)  = u_init
1377                v(:,j,i)  = v_init
1378             ENDDO
1379          ENDDO
1380!
1381!--       Mask topography
1382          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1383          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )
1384!
1385!--       Set initial horizontal velocities at the lowest computational grid
1386!--       levels to zero in order to avoid too small time steps caused by the
1387!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
1388!--       in the limiting formula!).
1389!--       Please note, in case land- or urban-surface model is used and a
1390!--       spinup is applied, masking the lowest computational level is not
1391!--       possible as MOST as well as energy-balance parametrizations will not
1392!--       work with zero wind velocity.
1393          IF ( ibc_uv_b /= 1  .AND.  .NOT.  spinup )  THEN
1394             DO  i = nxlg, nxrg
1395                DO  j = nysg, nyng
1396                   DO  k = nzb, nzt
1397                      u(k,j,i) = MERGE( u(k,j,i), 0.0_wp,                      &
1398                                        BTEST( wall_flags_0(k,j,i), 20 ) )
1399                      v(k,j,i) = MERGE( v(k,j,i), 0.0_wp,                      &
1400                                        BTEST( wall_flags_0(k,j,i), 21 ) )
1401                   ENDDO
1402                ENDDO
1403             ENDDO
1404          ENDIF
1405
1406          IF ( humidity )  THEN
1407             DO  i = nxlg, nxrg
1408                DO  j = nysg, nyng
1409                   q(:,j,i) = q_init
1410                ENDDO
1411             ENDDO
1412          ENDIF
1413         
1414          IF ( passive_scalar )  THEN
1415             DO  i = nxlg, nxrg
1416                DO  j = nysg, nyng
1417                   s(:,j,i) = s_init
1418                ENDDO
1419             ENDDO
1420          ENDIF
1421
1422!
1423!--       Compute initial temperature field and other constants used in case
1424!--       of a sloping surface
1425          IF ( sloping_surface )  CALL init_slope
1426!
1427!--       Initialize surface variables, e.g. friction velocity, momentum
1428!--       fluxes, etc.
1429          CALL init_surfaces
1430         
1431          IF ( debug_output )  CALL debug_message( 'initializing with constant profiles', 'end' )
1432
1433       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )                  &
1434       THEN
1435
1436          IF ( debug_output )  CALL debug_message( 'initializing by user', 'start' )
1437!
1438!--       Pre-initialize surface variables, i.e. setting start- and end-indices
1439!--       at each (j,i)-location. Please note, this does not supersede
1440!--       user-defined initialization of surface quantities.
1441          CALL init_surfaces
1442!
1443!--       Initialization will completely be done by the user
1444          CALL user_init_3d_model
1445
1446          IF ( debug_output )  CALL debug_message( 'initializing by user', 'end' )
1447
1448       ENDIF
1449
1450       IF ( debug_output )  CALL debug_message( 'initializing statistics, boundary conditions, etc.', 'start' )
1451
1452!
1453!--    Bottom boundary
1454       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
1455          u(nzb,:,:) = 0.0_wp
1456          v(nzb,:,:) = 0.0_wp
1457       ENDIF
1458
1459!
1460!--    Apply channel flow boundary condition
1461       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
1462          u(nzt+1,:,:) = 0.0_wp
1463          v(nzt+1,:,:) = 0.0_wp
1464       ENDIF
1465
1466!
1467!--    Calculate virtual potential temperature
1468       IF ( humidity )  vpt = pt * ( 1.0_wp + 0.61_wp * q )
1469
1470!
1471!--    Store initial profiles for output purposes etc.. Please note, in case of
1472!--    initialization of u, v, w, pt, and q via output data derived from larger
1473!--    scale models, data will not be horizontally homogeneous. Actually, a mean
1474!--    profile should be calculated before.   
1475       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
1476       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
1477       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
1478          hom(nzb,1,5,:) = 0.0_wp
1479          hom(nzb,1,6,:) = 0.0_wp
1480       ENDIF
1481       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1482
1483       IF ( humidity )  THEN
1484!
1485!--       Store initial profile of total water content, virtual potential
1486!--       temperature
1487          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
1488          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
1489!
1490!--       Store initial profile of mixing ratio and potential
1491!--       temperature
1492          IF ( bulk_cloud_model  .OR.  cloud_droplets ) THEN
1493             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1494             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1495          ENDIF
1496       ENDIF
1497
1498!
1499!--    Store initial scalar profile
1500       IF ( passive_scalar )  THEN
1501          hom(:,1,121,:) = SPREAD(  s(:,nys,nxl), 2, statistic_regions+1 )
1502       ENDIF
1503
1504!
1505!--    Initialize the random number generators (from numerical recipes)
1506       CALL random_function_ini
1507       
1508       IF ( random_generator == 'random-parallel' )  THEN
1509          CALL init_parallel_random_generator( nx, nys, nyn, nxl, nxr )
1510       ENDIF
1511!
1512!--    Set the reference state to be used in the buoyancy terms (for ocean runs
1513!--    the reference state will be set (overwritten) in init_ocean)
1514       IF ( use_single_reference_value )  THEN
1515          IF (  .NOT.  humidity )  THEN
1516             ref_state(:) = pt_reference
1517          ELSE
1518             ref_state(:) = vpt_reference
1519          ENDIF
1520       ELSE
1521          IF (  .NOT.  humidity )  THEN
1522             ref_state(:) = pt_init(:)
1523          ELSE
1524             ref_state(:) = vpt(:,nys,nxl)
1525          ENDIF
1526       ENDIF
1527
1528!
1529!--    For the moment, vertical velocity is zero
1530       w = 0.0_wp
1531
1532!
1533!--    Initialize array sums (must be defined in first call of pres)
1534       sums = 0.0_wp
1535
1536!
1537!--    In case of iterative solvers, p must get an initial value
1538       IF ( psolver(1:9) == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0_wp
1539!
1540!--    Impose vortex with vertical axis on the initial velocity profile
1541       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
1542          CALL init_rankine
1543       ENDIF
1544
1545!
1546!--    Impose temperature anomaly (advection test only) or warm air bubble
1547!--    close to surface
1548       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0  .OR.  &
1549            INDEX( initializing_actions, 'initialize_bubble' ) /= 0  )  THEN
1550          CALL init_pt_anomaly
1551       ENDIF
1552       
1553!
1554!--    If required, change the surface temperature at the start of the 3D run
1555       IF ( pt_surface_initial_change /= 0.0_wp )  THEN
1556          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
1557       ENDIF
1558
1559!
1560!--    If required, change the surface humidity/scalar at the start of the 3D
1561!--    run
1562       IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )              &
1563          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
1564         
1565       IF ( passive_scalar .AND.  s_surface_initial_change /= 0.0_wp )         &
1566          s(nzb,:,:) = s(nzb,:,:) + s_surface_initial_change
1567       
1568
1569!
1570!--    Initialize old and new time levels.
1571       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1572       pt_p = pt; u_p = u; v_p = v; w_p = w
1573
1574       IF ( humidity  )  THEN
1575          tq_m = 0.0_wp
1576          q_p = q
1577       ENDIF
1578       
1579       IF ( passive_scalar )  THEN
1580          ts_m = 0.0_wp
1581          s_p  = s
1582       ENDIF       
1583
1584       IF ( debug_output )  CALL debug_message( 'initializing statistics, boundary conditions, etc.', 'end' )
1585
1586    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.         &
1587             TRIM( initializing_actions ) == 'cyclic_fill' )                   &
1588    THEN
1589
1590       IF ( debug_output )  CALL debug_message( 'initializing in case of restart / cyclic_fill', 'start' )
1591!
1592!--    Initialize surface elements and its attributes, e.g. heat- and
1593!--    momentumfluxes, roughness, scaling parameters. As number of surface
1594!--    elements might be different between runs, e.g. in case of cyclic fill,
1595!--    and not all surface elements are read, surface elements need to be
1596!--    initialized before.
1597!--    Please note, in case of cyclic fill, surfaces should be initialized
1598!--    after restart data is read, else, individual settings of surface
1599!--    parameters will be overwritten from data of precursor run, hence,
1600!--    init_surfaces is called a second time after reading the restart data.
1601       CALL init_surfaces                       
1602!
1603!--    When reading data for cyclic fill of 3D prerun data files, read
1604!--    some of the global variables from the restart file which are required
1605!--    for initializing the inflow
1606       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1607
1608          DO  i = 0, io_blocks-1
1609             IF ( i == io_group )  THEN
1610                CALL rrd_read_parts_of_global
1611             ENDIF
1612#if defined( __parallel )
1613             CALL MPI_BARRIER( comm2d, ierr )
1614#endif
1615          ENDDO
1616
1617       ENDIF
1618
1619!
1620!--    Read processor specific binary data from restart file
1621       DO  i = 0, io_blocks-1
1622          IF ( i == io_group )  THEN
1623             CALL rrd_local
1624          ENDIF
1625#if defined( __parallel )
1626          CALL MPI_BARRIER( comm2d, ierr )
1627#endif
1628       ENDDO
1629!
1630!--    In case of cyclic fill, call init_surfaces a second time, so that
1631!--    surface properties such as heat fluxes are initialized as prescribed.
1632       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )                    &
1633          CALL init_surfaces
1634
1635!
1636!--    In case of complex terrain and cyclic fill method as initialization,
1637!--    shift initial data in the vertical direction for each point in the
1638!--    x-y-plane depending on local surface height
1639       IF ( complex_terrain  .AND.                                             &
1640            TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1641          DO  i = nxlg, nxrg
1642             DO  j = nysg, nyng
1643                nz_u_shift = get_topography_top_index_ji( j, i, 'u' )
1644                nz_v_shift = get_topography_top_index_ji( j, i, 'v' )
1645                nz_w_shift = get_topography_top_index_ji( j, i, 'w' )
1646                nz_s_shift = get_topography_top_index_ji( j, i, 's' )
1647
1648                u(nz_u_shift:nzt+1,j,i)  = u(0:nzt+1-nz_u_shift,j,i)               
1649
1650                v(nz_v_shift:nzt+1,j,i)  = v(0:nzt+1-nz_v_shift,j,i)
1651
1652                w(nz_w_shift:nzt+1,j,i)  = w(0:nzt+1-nz_w_shift,j,i)
1653
1654                p(nz_s_shift:nzt+1,j,i)  =  p(0:nzt+1-nz_s_shift,j,i)
1655                pt(nz_s_shift:nzt+1,j,i) = pt(0:nzt+1-nz_s_shift,j,i)
1656             ENDDO
1657          ENDDO
1658       ENDIF
1659
1660!
1661!--    Initialization of the turbulence recycling method
1662       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.               &
1663            turbulent_inflow )  THEN
1664!
1665!--       First store the profiles to be used at the inflow.
1666!--       These profiles are the (temporally) and horizontally averaged vertical
1667!--       profiles from the prerun. Alternatively, prescribed profiles
1668!--       for u,v-components can be used.
1669          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,1:num_mean_inflow_profiles) )
1670
1671          IF ( use_prescribed_profile_data )  THEN
1672             mean_inflow_profiles(:,1) = u_init            ! u
1673             mean_inflow_profiles(:,2) = v_init            ! v
1674          ELSE
1675             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
1676             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
1677          ENDIF
1678          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
1679          IF ( humidity )                                                      &
1680             mean_inflow_profiles(:,6) = hom_sum(:,41,0)   ! q
1681          IF ( passive_scalar )                                                &
1682             mean_inflow_profiles(:,7) = hom_sum(:,115,0)   ! s
1683!
1684!--       In case of complex terrain, determine vertical displacement at inflow
1685!--       boundary and adjust mean inflow profiles
1686          IF ( complex_terrain )  THEN
1687             IF ( nxlg <= 0 .AND. nxrg >= 0 .AND. nysg <= 0 .AND. nyng >= 0 )  THEN
1688                nz_u_shift_l = get_topography_top_index_ji( 0, 0, 'u' )
1689                nz_v_shift_l = get_topography_top_index_ji( 0, 0, 'v' )
1690                nz_w_shift_l = get_topography_top_index_ji( 0, 0, 'w' )
1691                nz_s_shift_l = get_topography_top_index_ji( 0, 0, 's' )
1692             ELSE
1693                nz_u_shift_l = 0
1694                nz_v_shift_l = 0
1695                nz_w_shift_l = 0
1696                nz_s_shift_l = 0
1697             ENDIF
1698
1699#if defined( __parallel )
1700             CALL MPI_ALLREDUCE(nz_u_shift_l, nz_u_shift, 1, MPI_INTEGER,      &
1701                                MPI_MAX, comm2d, ierr)
1702             CALL MPI_ALLREDUCE(nz_v_shift_l, nz_v_shift, 1, MPI_INTEGER,      &
1703                                MPI_MAX, comm2d, ierr)
1704             CALL MPI_ALLREDUCE(nz_w_shift_l, nz_w_shift, 1, MPI_INTEGER,      & 
1705                                MPI_MAX, comm2d, ierr)
1706             CALL MPI_ALLREDUCE(nz_s_shift_l, nz_s_shift, 1, MPI_INTEGER,      &
1707                                MPI_MAX, comm2d, ierr)
1708#else
1709             nz_u_shift = nz_u_shift_l
1710             nz_v_shift = nz_v_shift_l
1711             nz_w_shift = nz_w_shift_l
1712             nz_s_shift = nz_s_shift_l
1713#endif
1714
1715             mean_inflow_profiles(:,1) = 0.0_wp
1716             mean_inflow_profiles(nz_u_shift:nzt+1,1) = hom_sum(0:nzt+1-nz_u_shift,1,0)  ! u
1717
1718             mean_inflow_profiles(:,2) = 0.0_wp
1719             mean_inflow_profiles(nz_v_shift:nzt+1,2) = hom_sum(0:nzt+1-nz_v_shift,2,0)  ! v
1720
1721             mean_inflow_profiles(nz_s_shift:nzt+1,4) = hom_sum(0:nzt+1-nz_s_shift,4,0)  ! pt
1722
1723          ENDIF
1724
1725!
1726!--       If necessary, adjust the horizontal flow field to the prescribed
1727!--       profiles
1728          IF ( use_prescribed_profile_data )  THEN
1729             DO  i = nxlg, nxrg
1730                DO  j = nysg, nyng
1731                   DO  k = nzb, nzt+1
1732                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1733                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
1734                   ENDDO
1735                ENDDO
1736             ENDDO
1737          ENDIF
1738
1739!
1740!--       Use these mean profiles at the inflow (provided that Dirichlet
1741!--       conditions are used)
1742          IF ( bc_dirichlet_l )  THEN
1743             DO  j = nysg, nyng
1744                DO  k = nzb, nzt+1
1745                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1746                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
1747                   w(k,j,nxlg:-1)  = 0.0_wp
1748                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
1749                   IF ( humidity )                                             &
1750                      q(k,j,nxlg:-1)  = mean_inflow_profiles(k,6)
1751                   IF ( passive_scalar )                                       &
1752                      s(k,j,nxlg:-1)  = mean_inflow_profiles(k,7)                     
1753                ENDDO
1754             ENDDO
1755          ENDIF
1756
1757!
1758!--       Calculate the damping factors to be used at the inflow. For a
1759!--       turbulent inflow the turbulent fluctuations have to be limited
1760!--       vertically because otherwise the turbulent inflow layer will grow
1761!--       in time.
1762          IF ( inflow_damping_height == 9999999.9_wp )  THEN
1763!
1764!--          Default: use the inversion height calculated by the prerun; if
1765!--          this is zero, inflow_damping_height must be explicitly
1766!--          specified.
1767             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0_wp )  THEN
1768                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1769             ELSE
1770                WRITE( message_string, * ) 'inflow_damping_height must be ',   &
1771                     'explicitly specified because&the inversion height ',     &
1772                     'calculated by the prerun is zero.'
1773                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
1774             ENDIF
1775
1776          ENDIF
1777
1778          IF ( inflow_damping_width == 9999999.9_wp )  THEN
1779!
1780!--          Default for the transition range: one tenth of the undamped
1781!--          layer
1782             inflow_damping_width = 0.1_wp * inflow_damping_height
1783
1784          ENDIF
1785
1786          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
1787
1788          DO  k = nzb, nzt+1
1789
1790             IF ( zu(k) <= inflow_damping_height )  THEN
1791                inflow_damping_factor(k) = 1.0_wp
1792             ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) )  THEN
1793                inflow_damping_factor(k) = 1.0_wp -                            &
1794                                           ( zu(k) - inflow_damping_height ) / &
1795                                           inflow_damping_width
1796             ELSE
1797                inflow_damping_factor(k) = 0.0_wp
1798             ENDIF
1799
1800          ENDDO
1801
1802       ENDIF
1803
1804!
1805!--    Inside buildings set velocities back to zero
1806       IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND.                &
1807            topography /= 'flat' )  THEN
1808!
1809!--       Inside buildings set velocities back to zero.
1810!--       Other scalars (pt, q, s, p, sa, ...) are ignored at present,
1811!--       maybe revise later.
1812          DO  i = nxlg, nxrg
1813             DO  j = nysg, nyng
1814                DO  k = nzb, nzt
1815                   u(k,j,i)     = MERGE( u(k,j,i), 0.0_wp,                     &
1816                                         BTEST( wall_flags_0(k,j,i), 1 ) )
1817                   v(k,j,i)     = MERGE( v(k,j,i), 0.0_wp,                     &
1818                                         BTEST( wall_flags_0(k,j,i), 2 ) )
1819                   w(k,j,i)     = MERGE( w(k,j,i), 0.0_wp,                     &
1820                                         BTEST( wall_flags_0(k,j,i), 3 ) )
1821                ENDDO
1822             ENDDO
1823          ENDDO
1824
1825       ENDIF
1826
1827!
1828!--    Calculate initial temperature field and other constants used in case
1829!--    of a sloping surface
1830       IF ( sloping_surface )  CALL init_slope
1831
1832!
1833!--    Initialize new time levels (only done in order to set boundary values
1834!--    including ghost points)
1835       pt_p = pt; u_p = u; v_p = v; w_p = w
1836       IF ( humidity )  THEN
1837          q_p = q
1838       ENDIF
1839       IF ( passive_scalar )  s_p  = s
1840!
1841!--    Allthough tendency arrays are set in prognostic_equations, they have
1842!--    have to be predefined here because they are used (but multiplied with 0)
1843!--    there before they are set.
1844       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1845       IF ( humidity )  THEN
1846          tq_m = 0.0_wp
1847       ENDIF
1848       IF ( passive_scalar )  ts_m  = 0.0_wp
1849
1850       IF ( debug_output )  CALL debug_message( 'initializing in case of restart / cyclic_fill', 'end' )
1851
1852    ELSE
1853!
1854!--    Actually this part of the programm should not be reached
1855       message_string = 'unknown initializing problem'
1856       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
1857    ENDIF
1858
1859
1860    IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1861!
1862!--    Initialize old timelevels needed for radiation boundary conditions
1863       IF ( bc_radiation_l )  THEN
1864          u_m_l(:,:,:) = u(:,:,1:2)
1865          v_m_l(:,:,:) = v(:,:,0:1)
1866          w_m_l(:,:,:) = w(:,:,0:1)
1867       ENDIF
1868       IF ( bc_radiation_r )  THEN
1869          u_m_r(:,:,:) = u(:,:,nx-1:nx)
1870          v_m_r(:,:,:) = v(:,:,nx-1:nx)
1871          w_m_r(:,:,:) = w(:,:,nx-1:nx)
1872       ENDIF
1873       IF ( bc_radiation_s )  THEN
1874          u_m_s(:,:,:) = u(:,0:1,:)
1875          v_m_s(:,:,:) = v(:,1:2,:)
1876          w_m_s(:,:,:) = w(:,0:1,:)
1877       ENDIF
1878       IF ( bc_radiation_n )  THEN
1879          u_m_n(:,:,:) = u(:,ny-1:ny,:)
1880          v_m_n(:,:,:) = v(:,ny-1:ny,:)
1881          w_m_n(:,:,:) = w(:,ny-1:ny,:)
1882       ENDIF
1883       
1884    ENDIF
1885
1886!
1887!-- Calculate the initial volume flow at the right and north boundary
1888    IF ( conserve_volume_flow )  THEN
1889
1890       IF ( use_prescribed_profile_data )  THEN
1891
1892          volume_flow_initial_l = 0.0_wp
1893          volume_flow_area_l    = 0.0_wp
1894
1895          IF ( nxr == nx )  THEN
1896             DO  j = nys, nyn
1897                DO  k = nzb+1, nzt
1898                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1899                                              u_init(k) * dzw(k)               &
1900                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1901                                              BTEST( wall_flags_0(k,j,nxr), 1 )&
1902                                            )
1903
1904                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
1905                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1906                                              BTEST( wall_flags_0(k,j,nxr), 1 )&
1907                                            )
1908                ENDDO
1909             ENDDO
1910          ENDIF
1911         
1912          IF ( nyn == ny )  THEN
1913             DO  i = nxl, nxr
1914                DO  k = nzb+1, nzt
1915                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1916                                              v_init(k) * dzw(k)               &       
1917                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1918                                              BTEST( wall_flags_0(k,nyn,i), 2 )&
1919                                            )
1920                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
1921                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1922                                              BTEST( wall_flags_0(k,nyn,i), 2 )&
1923                                            )
1924                ENDDO
1925             ENDDO
1926          ENDIF
1927
1928#if defined( __parallel )
1929          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1930                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1931          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1932                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1933
1934#else
1935          volume_flow_initial = volume_flow_initial_l
1936          volume_flow_area    = volume_flow_area_l
1937#endif 
1938
1939       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1940
1941          volume_flow_initial_l = 0.0_wp
1942          volume_flow_area_l    = 0.0_wp
1943
1944          IF ( nxr == nx )  THEN
1945             DO  j = nys, nyn
1946                DO  k = nzb+1, nzt
1947                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1948                                              hom_sum(k,1,0) * dzw(k)          &
1949                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1950                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
1951                                            )
1952                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
1953                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1954                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
1955                                            )
1956                ENDDO
1957             ENDDO
1958          ENDIF
1959         
1960          IF ( nyn == ny )  THEN
1961             DO  i = nxl, nxr
1962                DO  k = nzb+1, nzt
1963                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1964                                              hom_sum(k,2,0) * dzw(k)          &       
1965                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1966                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
1967                                            )
1968                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
1969                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1970                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
1971                                            )
1972                ENDDO
1973             ENDDO
1974          ENDIF
1975
1976#if defined( __parallel )
1977          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1978                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1979          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1980                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1981
1982#else
1983          volume_flow_initial = volume_flow_initial_l
1984          volume_flow_area    = volume_flow_area_l
1985#endif 
1986
1987       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1988
1989          volume_flow_initial_l = 0.0_wp
1990          volume_flow_area_l    = 0.0_wp
1991
1992          IF ( nxr == nx )  THEN
1993             DO  j = nys, nyn
1994                DO  k = nzb+1, nzt
1995                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1996                                              u(k,j,nx) * dzw(k)               &
1997                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1998                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
1999                                            )
2000                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
2001                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2002                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
2003                                            )
2004                ENDDO
2005             ENDDO
2006          ENDIF
2007         
2008          IF ( nyn == ny )  THEN
2009             DO  i = nxl, nxr
2010                DO  k = nzb+1, nzt
2011                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
2012                                              v(k,ny,i) * dzw(k)               &       
2013                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2014                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
2015                                            )
2016                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
2017                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2018                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
2019                                            )
2020                ENDDO
2021             ENDDO
2022          ENDIF
2023
2024#if defined( __parallel )
2025          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
2026                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
2027          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
2028                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
2029
2030#else
2031          volume_flow_initial = volume_flow_initial_l
2032          volume_flow_area    = volume_flow_area_l
2033#endif 
2034
2035       ENDIF
2036
2037!
2038!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated
2039!--    from u|v_bulk instead
2040       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
2041          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
2042          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
2043       ENDIF
2044
2045    ENDIF
2046!
2047!-- Finally, if random_heatflux is set, disturb shf at horizontal
2048!-- surfaces. Actually, this should be done in surface_mod, where all other
2049!-- initializations of surface quantities are done. However, this
2050!-- would create a ring dependency, hence, it is done here. Maybe delete
2051!-- disturb_heatflux and tranfer the respective code directly into the
2052!-- initialization in surface_mod.         
2053    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
2054         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
2055 
2056       IF ( use_surface_fluxes  .AND.  constant_heatflux  .AND.                &
2057            random_heatflux )  THEN
2058          IF ( surf_def_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_def_h(0) )
2059          IF ( surf_lsm_h%ns    >= 1 )  CALL disturb_heatflux( surf_lsm_h    )
2060          IF ( surf_usm_h%ns    >= 1 )  CALL disturb_heatflux( surf_usm_h    )
2061       ENDIF
2062    ENDIF
2063
2064!
2065!-- Compute total sum of grid points and the mean surface level height for each
2066!-- statistic region. These are mainly used for horizontal averaging of
2067!-- turbulence statistics.
2068!-- ngp_2dh: number of grid points of a horizontal cross section through the
2069!--          respective statistic region
2070!-- ngp_3d:  number of grid points of the respective statistic region
2071    ngp_2dh_outer_l   = 0
2072    ngp_2dh_outer     = 0
2073    ngp_2dh_s_inner_l = 0
2074    ngp_2dh_s_inner   = 0
2075    ngp_2dh_l         = 0
2076    ngp_2dh           = 0
2077    ngp_3d_inner_l    = 0.0_wp
2078    ngp_3d_inner      = 0
2079    ngp_3d            = 0
2080    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
2081
2082    mean_surface_level_height   = 0.0_wp
2083    mean_surface_level_height_l = 0.0_wp
2084!
2085!-- To do: New concept for these non-topography grid points!
2086    DO  sr = 0, statistic_regions
2087       DO  i = nxl, nxr
2088          DO  j = nys, nyn
2089             IF ( rmask(j,i,sr) == 1.0_wp )  THEN
2090!
2091!--             All xy-grid points
2092                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
2093!
2094!--             Determine mean surface-level height. In case of downward-
2095!--             facing walls are present, more than one surface level exist.
2096!--             In this case, use the lowest surface-level height.
2097                IF ( surf_def_h(0)%start_index(j,i) <=                         &
2098                     surf_def_h(0)%end_index(j,i) )  THEN
2099                   m = surf_def_h(0)%start_index(j,i)
2100                   k = surf_def_h(0)%k(m)
2101                   mean_surface_level_height_l(sr) =                           &
2102                                       mean_surface_level_height_l(sr) + zw(k-1)
2103                ENDIF
2104                IF ( surf_lsm_h%start_index(j,i) <=                            &
2105                     surf_lsm_h%end_index(j,i) )  THEN
2106                   m = surf_lsm_h%start_index(j,i)
2107                   k = surf_lsm_h%k(m)
2108                   mean_surface_level_height_l(sr) =                           &
2109                                       mean_surface_level_height_l(sr) + zw(k-1)
2110                ENDIF
2111                IF ( surf_usm_h%start_index(j,i) <=                            &
2112                     surf_usm_h%end_index(j,i) )  THEN
2113                   m = surf_usm_h%start_index(j,i)
2114                   k = surf_usm_h%k(m)
2115                   mean_surface_level_height_l(sr) =                           &
2116                                       mean_surface_level_height_l(sr) + zw(k-1)
2117                ENDIF
2118
2119                k_surf = k - 1
2120
2121                DO  k = nzb, nzt+1
2122!
2123!--                xy-grid points above topography
2124                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr)     +         &
2125                                  MERGE( 1, 0, BTEST( wall_flags_0(k,j,i), 24 ) )
2126
2127                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) +         &
2128                                  MERGE( 1, 0, BTEST( wall_flags_0(k,j,i), 22 ) )
2129
2130                ENDDO
2131!
2132!--             All grid points of the total domain above topography
2133                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + ( nz - k_surf + 2 )
2134
2135
2136
2137             ENDIF
2138          ENDDO
2139       ENDDO
2140    ENDDO
2141
2142    sr = statistic_regions + 1
2143#if defined( __parallel )
2144    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2145    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,    &
2146                        comm2d, ierr )
2147    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2148    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr,   &
2149                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
2150    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2151    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),          &
2152                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
2153    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2154    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL,  &
2155                        MPI_SUM, comm2d, ierr )
2156    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
2157    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2158    CALL MPI_ALLREDUCE( mean_surface_level_height_l(0),                        &
2159                        mean_surface_level_height(0), sr, MPI_REAL,            &
2160                        MPI_SUM, comm2d, ierr )
2161    mean_surface_level_height = mean_surface_level_height / REAL( ngp_2dh )
2162#else
2163    ngp_2dh         = ngp_2dh_l
2164    ngp_2dh_outer   = ngp_2dh_outer_l
2165    ngp_2dh_s_inner = ngp_2dh_s_inner_l
2166    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
2167    mean_surface_level_height = mean_surface_level_height_l / REAL( ngp_2dh_l )
2168#endif
2169
2170    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * &
2171             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
2172
2173!
2174!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
2175!-- buoyancy, etc. A zero value will occur for cases where all grid points of
2176!-- the respective subdomain lie below the surface topography
2177    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   ) 
2178    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )),             &
2179                           ngp_3d_inner(:) )
2180    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) 
2181
2182    DEALLOCATE( mean_surface_level_height_l, ngp_2dh_l, ngp_2dh_outer_l,       &
2183                ngp_3d_inner_l, ngp_3d_inner_tmp )
2184!
2185!-- Initialize surface forcing corresponding to large-scale forcing. Therein,
2186!-- initialize heat-fluxes, etc. via datatype. Revise it later!
2187    IF ( large_scale_forcing .AND. lsf_surf )  THEN
2188       IF ( use_surface_fluxes  .AND.  constant_heatflux )  THEN
2189          CALL ls_forcing_surf ( simulated_time )
2190       ENDIF
2191    ENDIF
2192!
2193!-- Initializae 3D offline nesting in COSMO model and read data from
2194!-- external NetCDF file.
2195    IF ( nesting_offline )  CALL nesting_offl_init
2196!
2197!-- Initialize quantities for special advections schemes
2198    CALL init_advec
2199
2200!
2201!-- Impose random perturbation on the horizontal velocity field and then
2202!-- remove the divergences from the velocity field at the initial stage
2203    IF ( create_disturbances  .AND.  disturbance_energy_limit /= 0.0_wp  .AND. &
2204         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
2205         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
2206
2207       IF ( debug_output )  CALL debug_message( 'creating disturbances + applying pressure solver', 'start' )
2208!
2209!--    Needed for both disturb_field and pres
2210!$ACC DATA &
2211!$ACC CREATE(tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
2212!$ACC COPY(u(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
2213!$ACC COPY(v(nzb:nzt+1,nysg:nyng,nxlg:nxrg))
2214
2215       CALL disturb_field( 'u', tend, u )
2216       CALL disturb_field( 'v', tend, v )
2217
2218!$ACC DATA &
2219!$ACC CREATE(d(nzb+1:nzt,nys:nyn,nxl:nxr)) &
2220!$ACC COPY(w(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
2221!$ACC COPY(p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
2222!$ACC COPYIN(rho_air(nzb:nzt+1), rho_air_zw(nzb:nzt+1)) &
2223!$ACC COPYIN(ddzu(1:nzt+1), ddzw(1:nzt+1)) &
2224!$ACC COPYIN(wall_flags_0(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
2225!$ACC COPYIN(bc_h(0:1)) &
2226!$ACC COPYIN(bc_h(0)%i(1:bc_h(0)%ns)) &
2227!$ACC COPYIN(bc_h(0)%j(1:bc_h(0)%ns)) &
2228!$ACC COPYIN(bc_h(0)%k(1:bc_h(0)%ns)) &
2229!$ACC COPYIN(bc_h(1)%i(1:bc_h(1)%ns)) &
2230!$ACC COPYIN(bc_h(1)%j(1:bc_h(1)%ns)) &
2231!$ACC COPYIN(bc_h(1)%k(1:bc_h(1)%ns))
2232
2233       n_sor = nsor_ini
2234       CALL pres
2235       n_sor = nsor
2236
2237!$ACC END DATA
2238!$ACC END DATA
2239
2240       IF ( debug_output )  CALL debug_message( 'creating disturbances + applying pressure solver', 'end' )
2241
2242    ENDIF
2243
2244    IF ( .NOT. ocean_mode )  THEN
2245
2246       ALLOCATE( hyp(nzb:nzt+1) )
2247       ALLOCATE( d_exner(nzb:nzt+1) )
2248       ALLOCATE( exner(nzb:nzt+1) )
2249       ALLOCATE( hyrho(nzb:nzt+1) )
2250!
2251!--    Check temperature in case of too large domain height
2252       DO  k = nzb, nzt+1
2253          IF ( ( pt_surface * exner_function(surface_pressure * 100.0_wp) - g/c_p * zu(k) ) < 0.0_wp )  THEN
2254             WRITE( message_string, * )  'absolute temperature < 0.0 at zu(', k, &
2255                                         ') = ', zu(k)
2256             CALL message( 'init_3d_model', 'PA0142', 1, 2, 0, 6, 0 )
2257          ENDIF
2258       ENDDO
2259
2260!
2261!--    Calculate vertical profile of the hydrostatic pressure (hyp)
2262       hyp    = barometric_formula(zu, pt_surface * exner_function(surface_pressure * 100.0_wp), surface_pressure * 100.0_wp)
2263       d_exner = exner_function_invers(hyp)
2264       exner = 1.0_wp / exner_function_invers(hyp)
2265       hyrho  = ideal_gas_law_rho_pt(hyp, pt_init)
2266!
2267!--    Compute reference density
2268       rho_surface = ideal_gas_law_rho(surface_pressure * 100.0_wp, pt_surface * exner_function(surface_pressure * 100.0_wp))
2269
2270    ENDIF
2271
2272!
2273!-- If required, initialize particles
2274    IF ( agents_active )  CALL mas_init
2275!
2276!-- Initialization of synthetic turbulence generator
2277    IF ( use_syn_turb_gen )  CALL stg_init
2278!
2279!-- Initializing actions for all other modules
2280    CALL module_interface_init
2281!
2282!-- Initialize surface layer (done after LSM as roughness length are required
2283!-- for initialization
2284    IF ( constant_flux_layer )  CALL init_surface_layer_fluxes
2285!
2286!-- Initialize surface data output
2287    IF ( surface_output )  CALL surface_data_output_init
2288!
2289!-- Initialize the ws-scheme.   
2290    IF ( ws_scheme_sca .OR. ws_scheme_mom )  CALL ws_init
2291!
2292!-- Perform post-initializing checks for all other modules
2293    CALL module_interface_init_checks
2294
2295!
2296!-- Setting weighting factors for calculation of perturbation pressure
2297!-- and turbulent quantities from the RK substeps
2298    IF ( TRIM(timestep_scheme) == 'runge-kutta-3' )  THEN      ! for RK3-method
2299
2300       weight_substep(1) = 1._wp/6._wp
2301       weight_substep(2) = 3._wp/10._wp
2302       weight_substep(3) = 8._wp/15._wp
2303
2304       weight_pres(1)    = 1._wp/3._wp
2305       weight_pres(2)    = 5._wp/12._wp
2306       weight_pres(3)    = 1._wp/4._wp
2307
2308    ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' )  THEN  ! for RK2-method
2309
2310       weight_substep(1) = 1._wp/2._wp
2311       weight_substep(2) = 1._wp/2._wp
2312         
2313       weight_pres(1)    = 1._wp/2._wp
2314       weight_pres(2)    = 1._wp/2._wp       
2315
2316    ELSE                                     ! for Euler-method
2317
2318       weight_substep(1) = 1.0_wp     
2319       weight_pres(1)    = 1.0_wp                   
2320
2321    ENDIF
2322
2323!
2324!-- Initialize Rayleigh damping factors
2325    rdf    = 0.0_wp
2326    rdf_sc = 0.0_wp
2327    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
2328
2329       IF (  .NOT.  ocean_mode )  THEN
2330          DO  k = nzb+1, nzt
2331             IF ( zu(k) >= rayleigh_damping_height )  THEN
2332                rdf(k) = rayleigh_damping_factor *                             &
2333                      ( SIN( pi * 0.5_wp * ( zu(k) - rayleigh_damping_height ) &
2334                             / ( zu(nzt) - rayleigh_damping_height ) )         &
2335                      )**2
2336             ENDIF
2337          ENDDO
2338       ELSE
2339!
2340!--       In ocean mode, rayleigh damping is applied in the lower part of the
2341!--       model domain
2342          DO  k = nzt, nzb+1, -1
2343             IF ( zu(k) <= rayleigh_damping_height )  THEN
2344                rdf(k) = rayleigh_damping_factor *                             &
2345                      ( SIN( pi * 0.5_wp * ( rayleigh_damping_height - zu(k) ) &
2346                             / ( rayleigh_damping_height - zu(nzb+1) ) )       &
2347                      )**2
2348             ENDIF
2349          ENDDO
2350       ENDIF
2351
2352    ENDIF
2353    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
2354
2355!
2356!-- Initialize the starting level and the vertical smoothing factor used for
2357!-- the external pressure gradient
2358    dp_smooth_factor = 1.0_wp
2359    IF ( dp_external )  THEN
2360!
2361!--    Set the starting level dp_level_ind_b only if it has not been set before
2362!--    (e.g. in init_grid).
2363       IF ( dp_level_ind_b == 0 )  THEN
2364          ind_array = MINLOC( ABS( dp_level_b - zu ) )
2365          dp_level_ind_b = ind_array(1) - 1 + nzb 
2366                                        ! MINLOC uses lower array bound 1
2367       ENDIF
2368       IF ( dp_smooth )  THEN
2369          dp_smooth_factor(:dp_level_ind_b) = 0.0_wp
2370          DO  k = dp_level_ind_b+1, nzt
2371             dp_smooth_factor(k) = 0.5_wp * ( 1.0_wp + SIN( pi *               &
2372                        ( REAL( k - dp_level_ind_b, KIND=wp ) /                &
2373                          REAL( nzt - dp_level_ind_b, KIND=wp ) - 0.5_wp ) ) )
2374          ENDDO
2375       ENDIF
2376    ENDIF
2377
2378!
2379!-- Initialize damping zone for the potential temperature in case of
2380!-- non-cyclic lateral boundaries. The damping zone has the maximum value
2381!-- at the inflow boundary and decreases to zero at pt_damping_width.
2382    ptdf_x = 0.0_wp
2383    ptdf_y = 0.0_wp
2384    IF ( bc_lr_dirrad )  THEN
2385       DO  i = nxl, nxr
2386          IF ( ( i * dx ) < pt_damping_width )  THEN
2387             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5_wp *              &
2388                            REAL( pt_damping_width - i * dx, KIND=wp ) / (     &
2389                            REAL( pt_damping_width, KIND=wp ) ) ) )**2 
2390          ENDIF
2391       ENDDO
2392    ELSEIF ( bc_lr_raddir )  THEN
2393       DO  i = nxl, nxr
2394          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
2395             ptdf_x(i) = pt_damping_factor *                                   &
2396                         SIN( pi * 0.5_wp *                                    &
2397                                 ( ( i - nx ) * dx + pt_damping_width ) /      &
2398                                 REAL( pt_damping_width, KIND=wp ) )**2
2399          ENDIF
2400       ENDDO 
2401    ELSEIF ( bc_ns_dirrad )  THEN
2402       DO  j = nys, nyn
2403          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
2404             ptdf_y(j) = pt_damping_factor *                                   &
2405                         SIN( pi * 0.5_wp *                                    &
2406                                 ( ( j - ny ) * dy + pt_damping_width ) /      &
2407                                 REAL( pt_damping_width, KIND=wp ) )**2
2408          ENDIF
2409       ENDDO 
2410    ELSEIF ( bc_ns_raddir )  THEN
2411       DO  j = nys, nyn
2412          IF ( ( j * dy ) < pt_damping_width )  THEN
2413             ptdf_y(j) = pt_damping_factor *                                   &
2414                         SIN( pi * 0.5_wp *                                    &
2415                                ( pt_damping_width - j * dy ) /                &
2416                                REAL( pt_damping_width, KIND=wp ) )**2
2417          ENDIF
2418       ENDDO
2419    ENDIF
2420
2421!
2422!-- Input binary data file is not needed anymore. This line must be placed
2423!-- after call of user_init!
2424    CALL close_file( 13 )
2425!
2426!-- In case of nesting, put an barrier to assure that all parent and child
2427!-- domains finished initialization.
2428#if defined( __parallel )
2429    IF ( nested_run )  CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
2430#endif
2431
2432
2433    CALL location_message( 'model initialization', 'finished' )
2434
2435 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.