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

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

optimization of two-moments cloud physics

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