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

Last change on this file since 1322 was 1322, checked in by raasch, 10 years ago

REAL functions and a lot of REAL constants provided with KIND-attribute,
some small bugfixes

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