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

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

New:
---

GPU porting of pres, swap_timelevel. Adjustments of openACC directives.
Further porting of poisfft, which now runs completely on GPU without any
host/device data transfer for serial an parallel runs (but parallel runs
require data transfer before and after the MPI transpositions).
GPU-porting of tridiagonal solver:
tridiagonal routines split into extermal subroutines (instead using CONTAINS),
no distinction between parallel/non-parallel in poisfft and tridia any more,
tridia routines moved to end of file because of probable bug in PGI compiler
(otherwise "invalid device function" is indicated during runtime).
(cuda_fft_interfaces, fft_xy, flow_statistics, init_3d_model, palm, poisfft, pres, prognostic_equations, swap_timelevel, time_integration, transpose)
output of accelerator board information. (header)

optimization of tridia routines: constant elements and coefficients of tri are
stored in seperate arrays ddzuw and tric, last dimension of tri reduced from 5 to 2,
(init_grid, init_3d_model, modules, palm, poisfft)

poisfft_init is now called internally from poisfft,
(Makefile, Makefile_check, init_pegrid, poisfft, poisfft_hybrid)

CPU-time per grid point and timestep is output to CPU_MEASURES file
(cpu_statistics, modules, time_integration)

Changed:


resorting from/to array work changed, work now has 4 dimensions instead of 1 (transpose)
array diss allocated only if required (init_3d_model)

pressure boundary condition "Neumann+inhomo" removed from the code
(check_parameters, header, poisfft, poisfft_hybrid, pres)

Errors:


bugfix: dependency added for cuda_fft_interfaces (Makefile)
bugfix: CUDA fft plans adjusted for domain decomposition (before they always
used total domain) (fft_xy)

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