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

Last change on this file since 1329 was 1323, checked in by raasch, 11 years ago

last commit documented

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