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

Last change on this file since 2272 was 2270, checked in by maronga, 7 years ago

major revisions in land surface model

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