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

Last change on this file since 1986 was 1961, checked in by suehring, 8 years ago

last commit documented

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