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

Last change on this file since 1159 was 1159, checked in by fricke, 11 years ago

Bugfix: In case of non-cyclic lateral boundary conditions, Neumann boundary conditions for the velocity components at the outflow are in fact radiation boundary conditions using the maximum phase velocity that ensures numerical stability (CFL-condition).
Logical operator use_cmax is now used instead of bc_lr_dirneu/_neudir.

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