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

Last change on this file since 1357 was 1354, checked in by heinze, 11 years ago

last commit documented

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