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

Last change on this file since 1360 was 1360, checked in by hoffmann, 10 years ago

last commit documented

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