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

Last change on this file since 2298 was 2292, checked in by schwenkel, 7 years ago

implementation of new bulk microphysics scheme

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