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

Last change on this file since 1165 was 1160, checked in by fricke, 12 years ago

last commit documented

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