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

Last change on this file since 2284 was 2277, checked in by kanani, 8 years ago

code documentation and cleanup

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