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

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

two-moment cloud physics implemented

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