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

Last change on this file since 2263 was 2259, checked in by gronemeier, 7 years ago

Implemented synthetic turbulence generator

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