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

Last change on this file since 3240 was 3234, checked in by schwenkel, 6 years ago

The increase of dots_num in case of radiation or land surface model must be done before user_init is called

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