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

Last change on this file since 1341 was 1341, checked in by kanani, 10 years ago

last commit documented

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