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

Last change on this file since 1060 was 1054, checked in by hoffmann, 12 years ago

last commit documented

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