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

Last change on this file since 1097 was 1093, checked in by raasch, 12 years ago

last commit documented

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