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

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

last commit documented

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