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

Last change on this file since 2629 was 2618, checked in by suehring, 7 years ago

Reorder calls for initializing surface quantities and give example of user-defined initialization of surface quantities

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