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

Last change on this file since 1065 was 1065, checked in by hoffmann, 11 years ago

cloud physics: rain sedimentation and turbulence effects

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