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

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

last commit documented

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