source: palm/trunk/SOURCE/init_grid.f90 @ 1762

Last change on this file since 1762 was 1762, checked in by hellstea, 8 years ago

Introduction of nested domain system

  • Property svn:keywords set to Id
File size: 69.0 KB
Line 
1!> @file init_grid.f90
2!-------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2015 Leibniz Universitaet Hannover
17!-------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
21! Introduction of nested domain feature
22!
23! Former revisions:
24! -----------------
25! $Id: init_grid.f90 1762 2016-02-25 12:31:13Z hellstea $
26!
27! 1743 2016-01-13 10:23:51Z raasch
28! Bugfix for calculation of nzb_s_outer and nzb_u_outer at north boundary of
29! total domain
30!
31! 1691 2015-10-26 16:17:44Z maronga
32! Renamed prandtl_layer to constant_flux_layer.
33!
34! 1682 2015-10-07 23:56:08Z knoop
35! Code annotations made doxygen readable
36!
37! 1677 2015-10-02 13:25:23Z boeske
38! Bugfix: Ghost points are included in wall_flags_0 and wall_flags_00
39!
40! 1675 2015-10-02 08:28:59Z gronemeier
41! Bugfix: Definition of topography grid levels
42!
43! 1660 2015-09-21 08:15:16Z gronemeier
44! Bugfix: Definition of topography grid levels if vertical grid stretching
45!         starts below the maximum topography height.
46!
47! 1580 2015-04-10 13:43:49Z suehring
48! Bugfix: setting flags for 5th order scheme near buildings
49!
50! 1575 2015-03-27 09:56:27Z raasch
51! adjustments for psolver-queries
52!
53! 1557 2015-03-05 16:43:04Z suehring
54! Adjustment for monotoinic limiter
55!
56! 1418 2014-06-06 13:05:08Z fricke
57! Bugfix: Change if-condition for stretched grid in the ocean, with the old
58!          condition and a negative value for dz_stretch_level the condition
59!          was always true for the whole model domain
60!
61! 1409 2014-05-23 12:11:32Z suehring
62! Bugfix: set wall_flags_0 at inflow and outflow boundary also for i <= nxlu
63! j <= nysv
64!
65! 1353 2014-04-08 15:21:23Z heinze
66! REAL constants provided with KIND-attribute
67!
68! 1322 2014-03-20 16:38:49Z raasch
69! REAL constants defined as wp-kind
70!
71! 1320 2014-03-20 08:40:49Z raasch
72! ONLY-attribute added to USE-statements,
73! kind-parameters added to all INTEGER and REAL declaration statements,
74! kinds are defined in new module kinds,
75! revision history before 2012 removed,
76! comment fields (!:) to be used for variable explanations added to
77! all variable declaration statements
78!
79! 1221 2013-09-10 08:59:13Z raasch
80! wall_flags_00 introduced to hold bits 32-63,
81! additional 3D-flag arrays for replacing the 2D-index array nzb_s_inner in
82! loops optimized for openACC (pres + flow_statistics)
83!
84! 1092 2013-02-02 11:24:22Z raasch
85! unused variables removed
86!
87! 1069 2012-11-28 16:18:43Z maronga
88! bugfix: added coupling_char to TOPOGRAPHY_DATA to allow topography in the ocean
89!          model in case of coupled runs
90!
91! 1036 2012-10-22 13:43:42Z raasch
92! code put under GPL (PALM 3.9)
93!
94! 1015 2012-09-27 09:23:24Z raasch
95! lower index for calculating wall_flags_0 set to nzb_w_inner instead of
96! nzb_w_inner+1
97!
98! 996 2012-09-07 10:41:47Z raasch
99! little reformatting
100!
101! 978 2012-08-09 08:28:32Z fricke
102! Bugfix: nzb_max is set to nzt at non-cyclic lateral boundaries
103! Bugfix: Set wall_flags_0 for inflow boundary
104!
105! 927 2012-06-06 19:15:04Z raasch
106! Wall flags are not set for multigrid method in case of masking method
107!
108! 864 2012-03-27 15:10:33Z gryschka
109! In case of ocean and Dirichlet bottom bc for u and v dzu_mg and ddzu_pres
110! were not correctly defined for k=1.
111!
112! 861 2012-03-26 14:18:34Z suehring
113! Set wall_flags_0. The array is needed for degradation in ws-scheme near walls,
114! inflow and outflow boundaries as well as near the bottom and the top of the
115! model domain.!
116! Initialization of nzb_s_inner and nzb_w_inner.
117! gls has to be at least nbgp to do not exceed the array bounds of nzb_local
118! while setting wall_flags_0
119!
120! 843 2012-02-29 15:16:21Z gryschka
121! In case of ocean and dirichlet bc for u and v at the bottom
122! the first u-level ist defined at same height as the first w-level
123!
124! 818 2012-02-08 16:11:23Z maronga
125! Bugfix: topo_height is only required if topography is used. It is thus now
126! allocated in the topography branch
127!
128! 809 2012-01-30 13:32:58Z maronga
129! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
130!
131! 807 2012-01-25 11:53:51Z maronga
132! New cpp directive "__check" implemented which is used by check_namelist_files
133!
134! Revision 1.1  1997/08/11 06:17:45  raasch
135! Initial revision (Testversion)
136!
137!
138! Description:
139! ------------
140!> Creating grid depending constants
141!------------------------------------------------------------------------------!
142 SUBROUTINE init_grid
143 
144
145    USE arrays_3d,                                                             &
146        ONLY:  dd2zu, ddzu, ddzu_pres, ddzw, dzu, dzu_mg, dzw, dzw_mg, f1_mg,  &
147               f2_mg, f3_mg, l_grid, l_wall, zu, zw
148       
149    USE control_parameters,                                                    &
150        ONLY:  bc_lr, bc_ns, building_height, building_length_x,               &
151               building_length_y, building_wall_left, building_wall_south,     &
152               canyon_height, canyon_wall_left, canyon_wall_south,             &
153               canyon_width_x, canyon_width_y, constant_flux_layer,            &
154               coupling_char, dp_level_ind_b, dz, dz_max, dz_stretch_factor,   &
155               dz_stretch_level, dz_stretch_level_index, ibc_uv_b, io_blocks,  &
156               io_group, inflow_l, inflow_n, inflow_r, inflow_s,               &
157               masking_method, maximum_grid_level, message_string,             &
158               momentum_advec, nest_domain, nest_bound_l, nest_bound_n,        &
159               nest_bound_r, nest_bound_s, ocean, outflow_l, outflow_n,        &
160               outflow_r, outflow_s, psolver, scalar_advec, topography,        &
161               topography_grid_convention, use_surface_fluxes, use_top_fluxes, &
162               wall_adjustment_factor 
163       
164    USE grid_variables,                                                        &
165        ONLY:  ddx, ddx2, ddx2_mg, ddy, ddy2, ddy2_mg, dx, dx2, dy, dy2, fwxm, &
166               fwxp, fwym, fwyp, fxm, fxp, fym, fyp, wall_e_x, wall_e_y,       &
167               wall_u, wall_v, wall_w_x, wall_w_y, zu_s_inner, zw_w_inner
168       
169    USE indices,                                                               &
170        ONLY:  flags, nbgp, nx, nxl, nxlg, nxlu, nxl_mg, nxr, nxrg, nxr_mg,    &
171               ny, nyn, nyng, nyn_mg, nys, nysv, nys_mg, nysg, nz, nzb,        &
172               nzb_2d, nzb_diff, nzb_diff_s_inner, nzb_diff_s_outer,           &
173               nzb_diff_u, nzb_diff_v, nzb_max, nzb_s_inner, nzb_s_outer,      &
174               nzb_u_inner, nzb_u_outer, nzb_v_inner, nzb_v_outer,             &
175               nzb_w_inner, nzb_w_outer, nzt, nzt_diff, nzt_mg, rflags_invers, &
176               rflags_s_inner, wall_flags_0, wall_flags_00, wall_flags_1,      &
177               wall_flags_10, wall_flags_2, wall_flags_3,  wall_flags_4,       &
178               wall_flags_5, wall_flags_6, wall_flags_7, wall_flags_8,         &
179               wall_flags_9
180   
181    USE kinds
182   
183    USE pegrid
184
185    IMPLICIT NONE
186
187    INTEGER(iwp) ::  bh      !<
188    INTEGER(iwp) ::  blx     !<
189    INTEGER(iwp) ::  bly     !<
190    INTEGER(iwp) ::  bxl     !<
191    INTEGER(iwp) ::  bxr     !<
192    INTEGER(iwp) ::  byn     !<
193    INTEGER(iwp) ::  bys     !<
194    INTEGER(iwp) ::  ch      !<
195    INTEGER(iwp) ::  cwx     !<
196    INTEGER(iwp) ::  cwy     !<
197    INTEGER(iwp) ::  cxl     !<
198    INTEGER(iwp) ::  cxr     !<
199    INTEGER(iwp) ::  cyn     !<
200    INTEGER(iwp) ::  cys     !<
201    INTEGER(iwp) ::  gls     !<
202    INTEGER(iwp) ::  i       !<
203    INTEGER(iwp) ::  ii      !<
204    INTEGER(iwp) ::  inc     !<
205    INTEGER(iwp) ::  j       !<
206    INTEGER(iwp) ::  k       !<
207    INTEGER(iwp) ::  l       !<
208    INTEGER(iwp) ::  nxl_l   !<
209    INTEGER(iwp) ::  nxr_l   !<
210    INTEGER(iwp) ::  nyn_l   !<
211    INTEGER(iwp) ::  nys_l   !<
212    INTEGER(iwp) ::  nzb_si  !<
213    INTEGER(iwp) ::  nzt_l   !<
214    INTEGER(iwp) ::  vi      !<
215
216    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  vertical_influence  !<
217
218    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_nl  !<
219    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_nr  !<
220    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_sl  !<
221    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_sr  !<
222    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_l     !<
223    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_n     !<
224    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_r     !<
225    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_s     !<
226    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_local  !<
227    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_tmp    !<
228
229    LOGICAL :: flag_set = .FALSE.  !<
230
231    REAL(wp) ::  dx_l          !<
232    REAL(wp) ::  dy_l          !<
233    REAL(wp) ::  dz_stretched  !<
234
235    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  topo_height  !<
236
237   
238!
239!-- Calculation of horizontal array bounds including ghost layers
240    nxlg = nxl - nbgp
241    nxrg = nxr + nbgp
242    nysg = nys - nbgp
243    nyng = nyn + nbgp
244
245!
246!-- Allocate grid arrays
247    ALLOCATE( ddzu(1:nzt+1), ddzw(1:nzt+1), dd2zu(1:nzt), dzu(1:nzt+1),        &
248              dzw(1:nzt+1), l_grid(1:nzt), zu(nzb:nzt+1), zw(nzb:nzt+1) )
249
250!
251!-- Compute height of u-levels from constant grid length and dz stretch factors
252    IF ( dz == -1.0_wp )  THEN
253       message_string = 'missing dz'
254       CALL message( 'init_grid', 'PA0200', 1, 2, 0, 6, 0 ) 
255    ELSEIF ( dz <= 0.0_wp )  THEN
256       WRITE( message_string, * ) 'dz=',dz,' <= 0.0'
257       CALL message( 'init_grid', 'PA0201', 1, 2, 0, 6, 0 )
258    ENDIF
259
260!
261!-- Define the vertical grid levels
262    IF ( .NOT. ocean )  THEN
263!
264!--    Grid for atmosphere with surface at z=0 (k=0, w-grid).
265!--    The second u-level (k=1) corresponds to the top of the
266!--    Prandtl-layer.
267
268       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2 ) THEN
269          zu(0) = 0.0_wp
270      !    zu(0) = - dz * 0.5_wp
271       ELSE
272          zu(0) = - dz * 0.5_wp
273       ENDIF
274       zu(1) =   dz * 0.5_wp
275
276       dz_stretch_level_index = nzt+1
277       dz_stretched = dz
278       DO  k = 2, nzt+1
279          IF ( dz_stretch_level <= zu(k-1)  .AND.  dz_stretched < dz_max )  THEN
280             dz_stretched = dz_stretched * dz_stretch_factor
281             dz_stretched = MIN( dz_stretched, dz_max )
282             IF ( dz_stretch_level_index == nzt+1 ) dz_stretch_level_index = k-1
283          ENDIF
284          zu(k) = zu(k-1) + dz_stretched
285       ENDDO
286
287!
288!--    Compute the w-levels. They are always staggered half-way between the
289!--    corresponding u-levels. In case of dirichlet bc for u and v at the
290!--    ground the first u- and w-level (k=0) are defined at same height (z=0).
291!--    The top w-level is extrapolated linearly.
292       zw(0) = 0.0_wp
293       DO  k = 1, nzt
294          zw(k) = ( zu(k) + zu(k+1) ) * 0.5_wp
295       ENDDO
296       zw(nzt+1) = zw(nzt) + 2.0_wp * ( zu(nzt+1) - zw(nzt) )
297
298    ELSE
299!
300!--    Grid for ocean with free water surface is at k=nzt (w-grid).
301!--    In case of neumann bc at the ground the first first u-level (k=0) lies
302!--    below the first w-level (k=0). In case of dirichlet bc the first u- and
303!--    w-level are defined at same height, but staggered from the second level.
304!--    The second u-level (k=1) corresponds to the top of the Prandtl-layer.
305       zu(nzt+1) =   dz * 0.5_wp
306       zu(nzt)   = - dz * 0.5_wp
307
308       dz_stretch_level_index = 0
309       dz_stretched = dz
310       DO  k = nzt-1, 0, -1
311!
312!--       The default value of dz_stretch_level is positive, thus the first
313!--       condition is always true. Hence, the second condition is necessary.
314          IF ( dz_stretch_level >= zu(k+1)  .AND.  dz_stretch_level <= 0.0  &
315               .AND.  dz_stretched < dz_max )  THEN
316             dz_stretched = dz_stretched * dz_stretch_factor
317             dz_stretched = MIN( dz_stretched, dz_max )
318             IF ( dz_stretch_level_index == 0 ) dz_stretch_level_index = k+1
319          ENDIF
320          zu(k) = zu(k+1) - dz_stretched
321       ENDDO
322
323!
324!--    Compute the w-levels. They are always staggered half-way between the
325!--    corresponding u-levels, except in case of dirichlet bc for u and v
326!--    at the ground. In this case the first u- and w-level are defined at
327!--    same height. The top w-level (nzt+1) is not used but set for
328!--    consistency, since w and all scalar variables are defined up tp nzt+1.
329       zw(nzt+1) = dz
330       zw(nzt)   = 0.0_wp
331       DO  k = 0, nzt
332          zw(k) = ( zu(k) + zu(k+1) ) * 0.5_wp
333       ENDDO
334
335!
336!--    In case of dirichlet bc for u and v the first u- and w-level are defined
337!--    at same height.
338       IF ( ibc_uv_b == 0 ) THEN
339          zu(0) = zw(0)
340       ENDIF
341
342    ENDIF
343
344!
345!-- Compute grid lengths.
346    DO  k = 1, nzt+1
347       dzu(k)  = zu(k) - zu(k-1)
348       ddzu(k) = 1.0_wp / dzu(k)
349       dzw(k)  = zw(k) - zw(k-1)
350       ddzw(k) = 1.0_wp / dzw(k)
351    ENDDO
352
353    DO  k = 1, nzt
354       dd2zu(k) = 1.0_wp / ( dzu(k) + dzu(k+1) )
355    ENDDO
356   
357!   
358!-- The FFT- SOR-pressure solvers assume grid spacings of a staggered grid
359!-- everywhere. For the actual grid, the grid spacing at the lowest level
360!-- is only dz/2, but should be dz. Therefore, an additional array
361!-- containing with appropriate grid information is created for these
362!-- solvers.
363    IF ( psolver(1:9) /= 'multigrid' )  THEN
364       ALLOCATE( ddzu_pres(1:nzt+1) )
365       ddzu_pres = ddzu
366       ddzu_pres(1) = ddzu_pres(2)  ! change for lowest level
367    ENDIF   
368
369!
370!-- In case of multigrid method, compute grid lengths and grid factors for the
371!-- grid levels
372    IF ( psolver(1:9) == 'multigrid' )  THEN
373
374       ALLOCATE( ddx2_mg(maximum_grid_level), ddy2_mg(maximum_grid_level), &
375                 dzu_mg(nzb+1:nzt+1,maximum_grid_level),                   &
376                 dzw_mg(nzb+1:nzt+1,maximum_grid_level),                   &
377                 f1_mg(nzb+1:nzt,maximum_grid_level),                      &
378                 f2_mg(nzb+1:nzt,maximum_grid_level),                      &
379                 f3_mg(nzb+1:nzt,maximum_grid_level) )
380
381       dzu_mg(:,maximum_grid_level) = dzu
382!       
383!--    Next line to ensure an equally spaced grid.
384       dzu_mg(1,maximum_grid_level) = dzu(2)
385
386       dzw_mg(:,maximum_grid_level) = dzw
387       nzt_l = nzt
388       DO  l = maximum_grid_level-1, 1, -1
389           dzu_mg(nzb+1,l) = 2.0_wp * dzu_mg(nzb+1,l+1)
390           dzw_mg(nzb+1,l) = 2.0_wp * dzw_mg(nzb+1,l+1)
391           nzt_l = nzt_l / 2
392           DO  k = 2, nzt_l+1
393              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
394              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
395           ENDDO
396       ENDDO
397
398       nzt_l = nzt
399       dx_l  = dx
400       dy_l  = dy
401       DO  l = maximum_grid_level, 1, -1
402          ddx2_mg(l) = 1.0_wp / dx_l**2
403          ddy2_mg(l) = 1.0_wp / dy_l**2
404          DO  k = nzb+1, nzt_l
405             f2_mg(k,l) = 1.0_wp / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
406             f3_mg(k,l) = 1.0_wp / ( dzu_mg(k,l)   * dzw_mg(k,l) )
407             f1_mg(k,l) = 2.0_wp * ( ddx2_mg(l) + ddy2_mg(l) ) + &
408                          f2_mg(k,l) + f3_mg(k,l)
409          ENDDO
410          nzt_l = nzt_l / 2
411          dx_l  = dx_l * 2.0_wp
412          dy_l  = dy_l * 2.0_wp
413       ENDDO
414
415    ENDIF
416
417!
418!-- Compute the reciprocal values of the horizontal grid lengths.
419    ddx = 1.0_wp / dx
420    ddy = 1.0_wp / dy
421    dx2 = dx * dx
422    dy2 = dy * dy
423    ddx2 = 1.0_wp / dx2
424    ddy2 = 1.0_wp / dy2
425
426!
427!-- Compute the grid-dependent mixing length.
428    DO  k = 1, nzt
429       l_grid(k)  = ( dx * dy * dzw(k) )**0.33333333333333_wp
430    ENDDO
431
432!
433!-- Allocate outer and inner index arrays for topography and set
434!-- defaults.
435!-- nzb_local has to contain additional layers of ghost points for calculating
436!-- the flag arrays needed for the multigrid method
437    gls = 2**( maximum_grid_level )
438    IF ( gls < nbgp )  gls = nbgp
439
440    ALLOCATE( corner_nl(nys:nyn,nxl:nxr), corner_nr(nys:nyn,nxl:nxr),       &
441              corner_sl(nys:nyn,nxl:nxr), corner_sr(nys:nyn,nxl:nxr),       &
442              nzb_local(-gls:ny+gls,-gls:nx+gls),                                   &
443              nzb_tmp(-nbgp:ny+nbgp,-nbgp:nx+nbgp),                         &
444              wall_l(nys:nyn,nxl:nxr), wall_n(nys:nyn,nxl:nxr),             &
445              wall_r(nys:nyn,nxl:nxr), wall_s(nys:nyn,nxl:nxr) )
446    ALLOCATE( fwxm(nysg:nyng,nxlg:nxrg), fwxp(nysg:nyng,nxlg:nxrg),         &
447              fwym(nysg:nyng,nxlg:nxrg), fwyp(nysg:nyng,nxlg:nxrg),         &
448              fxm(nysg:nyng,nxlg:nxrg), fxp(nysg:nyng,nxlg:nxrg),           &
449              fym(nysg:nyng,nxlg:nxrg), fyp(nysg:nyng,nxlg:nxrg),           &
450              nzb_s_inner(nysg:nyng,nxlg:nxrg),                             &
451              nzb_s_outer(nysg:nyng,nxlg:nxrg),                             &
452              nzb_u_inner(nysg:nyng,nxlg:nxrg),                             &
453              nzb_u_outer(nysg:nyng,nxlg:nxrg),                             &
454              nzb_v_inner(nysg:nyng,nxlg:nxrg),                             &
455              nzb_v_outer(nysg:nyng,nxlg:nxrg),                             &
456              nzb_w_inner(nysg:nyng,nxlg:nxrg),                             &
457              nzb_w_outer(nysg:nyng,nxlg:nxrg),                             &
458              nzb_diff_s_inner(nysg:nyng,nxlg:nxrg),                        &
459              nzb_diff_s_outer(nysg:nyng,nxlg:nxrg),                        &
460              nzb_diff_u(nysg:nyng,nxlg:nxrg),                              &
461              nzb_diff_v(nysg:nyng,nxlg:nxrg),                              &
462              nzb_2d(nysg:nyng,nxlg:nxrg),                                  &
463              rflags_s_inner(nzb:nzt+2,nysg:nyng,nxlg:nxrg),                &
464              rflags_invers(nysg:nyng,nxlg:nxrg,nzb:nzt+2),                 &
465              wall_e_x(nysg:nyng,nxlg:nxrg),                                &
466              wall_e_y(nysg:nyng,nxlg:nxrg),                                &
467              wall_u(nysg:nyng,nxlg:nxrg),                                  &
468              wall_v(nysg:nyng,nxlg:nxrg),                                  &
469              wall_w_x(nysg:nyng,nxlg:nxrg),                                &
470              wall_w_y(nysg:nyng,nxlg:nxrg) )
471
472
473
474    ALLOCATE( l_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
475
476
477    nzb_s_inner = nzb;  nzb_s_outer = nzb
478    nzb_u_inner = nzb;  nzb_u_outer = nzb
479    nzb_v_inner = nzb;  nzb_v_outer = nzb
480    nzb_w_inner = nzb;  nzb_w_outer = nzb
481
482    rflags_s_inner = 1.0_wp
483    rflags_invers  = 1.0_wp
484
485!
486!-- Define vertical gridpoint from (or to) which on the usual finite difference
487!-- form (which does not use surface fluxes) is applied
488    IF ( constant_flux_layer  .OR.  use_surface_fluxes )  THEN
489       nzb_diff = nzb + 2
490    ELSE
491       nzb_diff = nzb + 1
492    ENDIF
493    IF ( use_top_fluxes )  THEN
494       nzt_diff = nzt - 1
495    ELSE
496       nzt_diff = nzt
497    ENDIF
498
499    nzb_diff_s_inner = nzb_diff;  nzb_diff_s_outer = nzb_diff
500    nzb_diff_u = nzb_diff;  nzb_diff_v = nzb_diff
501
502    wall_e_x = 0.0_wp;  wall_e_y = 0.0_wp;  wall_u = 0.0_wp;  wall_v = 0.0_wp
503    wall_w_x = 0.0_wp;  wall_w_y = 0.0_wp
504    fwxp = 1.0_wp;  fwxm = 1.0_wp;  fwyp = 1.0_wp;  fwym = 1.0_wp
505    fxp  = 1.0_wp;  fxm  = 1.0_wp;  fyp  = 1.0_wp;  fym  = 1.0_wp
506
507!
508!-- Initialize near-wall mixing length l_wall only in the vertical direction
509!-- for the moment,
510!-- multiplication with wall_adjustment_factor near the end of this routine
511    l_wall(nzb,:,:)   = l_grid(1)
512    DO  k = nzb+1, nzt
513       l_wall(k,:,:)  = l_grid(k)
514    ENDDO
515    l_wall(nzt+1,:,:) = l_grid(nzt)
516
517    ALLOCATE ( vertical_influence(nzb:nzt) )
518    DO  k = 1, nzt
519       vertical_influence(k) = MIN ( INT( l_grid(k) / &
520                     ( wall_adjustment_factor * dzw(k) ) + 0.5_wp ), nzt - k )
521    ENDDO
522
523    DO  k = 1, MAXVAL( nzb_s_inner )
524       IF ( l_grid(k) > 1.5_wp * dx * wall_adjustment_factor .OR.  &
525            l_grid(k) > 1.5_wp * dy * wall_adjustment_factor )  THEN
526          WRITE( message_string, * ) 'grid anisotropy exceeds ', &
527                                     'threshold given by only local', &
528                                     ' &horizontal reduction of near_wall ', &
529                                     'mixing length l_wall', &
530                                     ' &starting from height level k = ', k, '.'
531          CALL message( 'init_grid', 'PA0202', 0, 1, 0, 6, 0 )
532          EXIT
533       ENDIF
534    ENDDO
535    vertical_influence(0) = vertical_influence(1)
536
537    DO  i = nxlg, nxrg
538       DO  j = nysg, nyng
539          DO  k = nzb_s_inner(j,i) + 1, &
540                  nzb_s_inner(j,i) + vertical_influence(nzb_s_inner(j,i))
541             l_wall(k,j,i) = zu(k) - zw(nzb_s_inner(j,i))
542          ENDDO
543       ENDDO
544    ENDDO
545
546!
547!-- Set outer and inner index arrays for non-flat topography.
548!-- Here consistency checks concerning domain size and periodicity are
549!-- necessary.
550!-- Within this SELECT CASE structure only nzb_local is initialized
551!-- individually depending on the chosen topography type, all other index
552!-- arrays are initialized further below.
553    SELECT CASE ( TRIM( topography ) )
554
555       CASE ( 'flat' )
556!
557!--       nzb_local is required for the multigrid solver
558          nzb_local = 0
559
560       CASE ( 'single_building' )
561!
562!--       Single rectangular building, by default centered in the middle of the
563!--       total domain
564          blx = NINT( building_length_x / dx )
565          bly = NINT( building_length_y / dy )
566          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
567          IF ( ABS( zw(bh  ) - building_height ) == &
568               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
569
570          IF ( building_wall_left == 9999999.9_wp )  THEN
571             building_wall_left = ( nx + 1 - blx ) / 2 * dx
572          ENDIF
573          bxl = NINT( building_wall_left / dx )
574          bxr = bxl + blx
575
576          IF ( building_wall_south == 9999999.9_wp )  THEN
577             building_wall_south = ( ny + 1 - bly ) / 2 * dy
578          ENDIF
579          bys = NINT( building_wall_south / dy )
580          byn = bys + bly
581
582!
583!--       Building size has to meet some requirements
584          IF ( ( bxl < 1 ) .OR. ( bxr > nx-1 ) .OR. ( bxr < bxl+3 ) .OR.  &
585               ( bys < 1 ) .OR. ( byn > ny-1 ) .OR. ( byn < bys+3 ) )  THEN
586             WRITE( message_string, * ) 'inconsistent building parameters:',   &
587                                      '& bxl=', bxl, 'bxr=', bxr, 'bys=', bys, &
588                                      'byn=', byn, 'nx=', nx, 'ny=', ny
589             CALL message( 'init_grid', 'PA0203', 1, 2, 0, 6, 0 )
590          ENDIF
591
592!
593!--       Define the building.
594          nzb_local = 0
595          nzb_local(bys:byn,bxl:bxr) = bh
596
597       CASE ( 'single_street_canyon' )
598!
599!--       Single quasi-2D street canyon of infinite length in x or y direction.
600!--       The canyon is centered in the other direction by default.
601          IF ( canyon_width_x /= 9999999.9_wp )  THEN
602!
603!--          Street canyon in y direction
604             cwx = NINT( canyon_width_x / dx )
605             IF ( canyon_wall_left == 9999999.9_wp )  THEN
606                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
607             ENDIF
608             cxl = NINT( canyon_wall_left / dx )
609             cxr = cxl + cwx
610
611          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
612!
613!--          Street canyon in x direction
614             cwy = NINT( canyon_width_y / dy )
615             IF ( canyon_wall_south == 9999999.9_wp )  THEN
616                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
617             ENDIF
618             cys = NINT( canyon_wall_south / dy )
619             cyn = cys + cwy
620
621          ELSE
622             
623             message_string = 'no street canyon width given'
624             CALL message( 'init_grid', 'PA0204', 1, 2, 0, 6, 0 )
625 
626          ENDIF
627
628          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
629          IF ( ABS( zw(ch  ) - canyon_height ) == &
630               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
631
632          dp_level_ind_b = ch
633!
634!--       Street canyon size has to meet some requirements
635          IF ( canyon_width_x /= 9999999.9_wp )  THEN
636             IF ( ( cxl < 1 ) .OR. ( cxr > nx-1 ) .OR. ( cwx < 3 ) .OR.        &
637               ( ch < 3 ) )  THEN
638                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
639                                           '&cxl=', cxl, 'cxr=', cxr,          &
640                                           'cwx=', cwx,                        &
641                                           'ch=', ch, 'nx=', nx, 'ny=', ny
642                CALL message( 'init_grid', 'PA0205', 1, 2, 0, 6, 0 ) 
643             ENDIF
644          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
645             IF ( ( cys < 1 ) .OR. ( cyn > ny-1 ) .OR. ( cwy < 3 ) .OR.        &
646               ( ch < 3 ) )  THEN
647                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
648                                           '&cys=', cys, 'cyn=', cyn,          &
649                                           'cwy=', cwy,                        &
650                                           'ch=', ch, 'nx=', nx, 'ny=', ny
651                CALL message( 'init_grid', 'PA0206', 1, 2, 0, 6, 0 ) 
652             ENDIF
653          ENDIF
654          IF ( canyon_width_x /= 9999999.9_wp .AND.                            &                 
655               canyon_width_y /= 9999999.9_wp )  THEN
656             message_string = 'inconsistent canyon parameters:' //             &   
657                              '&street canyon can only be oriented' //         &
658                              '&either in x- or in y-direction'
659             CALL message( 'init_grid', 'PA0207', 1, 2, 0, 6, 0 )
660          ENDIF
661
662          nzb_local = ch
663          IF ( canyon_width_x /= 9999999.9_wp )  THEN
664             nzb_local(:,cxl+1:cxr-1) = 0
665          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
666             nzb_local(cys+1:cyn-1,:) = 0
667          ENDIF
668
669       CASE ( 'read_from_file' )
670
671          ALLOCATE ( topo_height(0:ny,0:nx) )
672
673          DO  ii = 0, io_blocks-1
674             IF ( ii == io_group )  THEN
675
676!
677!--             Arbitrary irregular topography data in PALM format (exactly
678!--             matching the grid size and total domain size)
679                OPEN( 90, FILE='TOPOGRAPHY_DATA'//coupling_char, STATUS='OLD', &
680                      FORM='FORMATTED', ERR=10 )
681                DO  j = ny, 0, -1
682                   READ( 90, *, ERR=11, END=11 )  ( topo_height(j,i), i = 0,nx )
683                ENDDO
684
685                GOTO 12
686         
687 10             message_string = 'file TOPOGRAPHY'//coupling_char//' does not exist'
688                CALL message( 'init_grid', 'PA0208', 1, 2, 0, 6, 0 )
689
690 11             message_string = 'errors in file TOPOGRAPHY_DATA'//coupling_char
691                CALL message( 'init_grid', 'PA0209', 1, 2, 0, 6, 0 )
692
693 12             CLOSE( 90 )
694
695             ENDIF
696#if defined( __parallel ) && ! defined ( __check )
697             CALL MPI_BARRIER( comm2d, ierr )
698#endif
699          ENDDO
700
701!
702!--       Calculate the index height of the topography
703          DO  i = 0, nx
704             DO  j = 0, ny
705                nzb_local(j,i) = MINLOC( ABS( zw - topo_height(j,i) ), 1 ) - 1
706                IF ( ABS( zw(nzb_local(j,i)  ) - topo_height(j,i) ) == &
707                     ABS( zw(nzb_local(j,i)+1) - topo_height(j,i) )    )  &
708                   nzb_local(j,i) = nzb_local(j,i) + 1
709             ENDDO
710          ENDDO
711
712          DEALLOCATE ( topo_height )
713!
714!--       Add cyclic boundaries (additional layers are for calculating
715!--       flag arrays needed for the multigrid sover)
716          nzb_local(-gls:-1,0:nx)     = nzb_local(ny-gls+1:ny,0:nx)
717          nzb_local(ny+1:ny+gls,0:nx) = nzb_local(0:gls-1,0:nx)
718          nzb_local(:,-gls:-1)        = nzb_local(:,nx-gls+1:nx)
719          nzb_local(:,nx+1:nx+gls)    = nzb_local(:,0:gls-1)
720
721       CASE DEFAULT
722!
723!--       The DEFAULT case is reached either if the parameter topography
724!--       contains a wrong character string or if the user has defined a special
725!--       case in the user interface. There, the subroutine user_init_grid
726!--       checks which of these two conditions applies.
727          CALL user_init_grid( gls, nzb_local )
728
729    END SELECT
730!
731!-- Determine the maximum level of topography. Furthermore it is used for
732!-- steering the degradation of order of the applied advection scheme.
733!-- In case of non-cyclic lateral boundaries, the order of the advection
734!-- scheme have to be reduced up to nzt (required at the lateral boundaries).
735    nzb_max = MAXVAL( nzb_local ) + 1
736    IF ( inflow_l .OR. outflow_l .OR. inflow_r .OR. outflow_r .OR.             &
737         inflow_n .OR. outflow_n .OR. inflow_s .OR. outflow_s .OR.             &
738         nest_domain )                                                         &
739    THEN
740       nzb_max = nzt
741    ENDIF
742
743!
744!-- Consistency checks and index array initialization are only required for
745!-- non-flat topography, also the initialization of topography height arrays
746!-- zu_s_inner and zw_w_inner
747    IF ( TRIM( topography ) /= 'flat' )  THEN
748
749!
750!--    Consistency checks
751       IF ( MINVAL( nzb_local ) < 0  .OR.  MAXVAL( nzb_local ) > nz + 1 )  THEN
752          WRITE( message_string, * ) 'nzb_local values are outside the',       &
753                                'model domain',                                &
754                                '&MINVAL( nzb_local ) = ', MINVAL(nzb_local),  &
755                                '&MAXVAL( nzb_local ) = ', MAXVAL(nzb_local)
756          CALL message( 'init_grid', 'PA0210', 1, 2, 0, 6, 0 )
757       ENDIF
758
759       IF ( bc_lr == 'cyclic' )  THEN
760          IF ( ANY( nzb_local(:,-1) /= nzb_local(:,nx)   )  .OR.               &
761               ANY( nzb_local(:,0)  /= nzb_local(:,nx+1) ) )  THEN
762             message_string = 'nzb_local does not fulfill cyclic' //           &
763                              ' boundary condition in x-direction'
764             CALL message( 'init_grid', 'PA0211', 1, 2, 0, 6, 0 )
765          ENDIF
766       ENDIF
767       IF ( bc_ns == 'cyclic' )  THEN
768          IF ( ANY( nzb_local(-1,:) /= nzb_local(ny,:)   )  .OR.               &
769               ANY( nzb_local(0,:)  /= nzb_local(ny+1,:) ) )  THEN
770             message_string = 'nzb_local does not fulfill cyclic' //           &
771                              ' boundary condition in y-direction'
772             CALL message( 'init_grid', 'PA0212', 1, 2, 0, 6, 0 )
773          ENDIF
774       ENDIF
775
776       IF ( topography_grid_convention == 'cell_edge' )  THEN
777!
778!--       The array nzb_local as defined using the 'cell_edge' convention
779!--       describes the actual total size of topography which is defined at the
780!--       cell edges where u=0 on the topography walls in x-direction and v=0
781!--       on the topography walls in y-direction. However, PALM uses individual
782!--       arrays nzb_u|v|w|s_inner|outer that are based on nzb_s_inner.
783!--       Therefore, the extent of topography in nzb_local is now reduced by
784!--       1dx at the E topography walls and by 1dy at the N topography walls
785!--       to form the basis for nzb_s_inner.
786          DO  j = -gls, ny + gls
787             DO  i = -gls, nx
788                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j,i+1) )
789             ENDDO
790          ENDDO
791!--       apply cyclic boundary conditions in x-direction
792!(ist das erforderlich? Ursache von Seung Bus Fehler?)
793          nzb_local(:,nx+1:nx+gls) = nzb_local(:,0:gls-1)
794          DO  i = -gls, nx + gls
795             DO  j = -gls, ny
796                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j+1,i) )
797             ENDDO
798          ENDDO
799!--       apply cyclic boundary conditions in y-direction
800!(ist das erforderlich? Ursache von Seung Bus Fehler?)
801          nzb_local(ny+1:ny+gls,:) = nzb_local(0:gls-1,:)
802       ENDIF
803
804!
805!--    Initialize index arrays nzb_s_inner and nzb_w_inner
806       nzb_s_inner = nzb_local(nysg:nyng,nxlg:nxrg)
807       nzb_w_inner = nzb_local(nysg:nyng,nxlg:nxrg)
808
809!
810!--    Initialize remaining index arrays:
811!--    first pre-initialize them with nzb_s_inner...
812       nzb_u_inner = nzb_s_inner
813       nzb_u_outer = nzb_s_inner
814       nzb_v_inner = nzb_s_inner
815       nzb_v_outer = nzb_s_inner
816       nzb_w_outer = nzb_s_inner
817       nzb_s_outer = nzb_s_inner
818
819!
820!--    ...then extend pre-initialized arrays in their according directions
821!--    based on nzb_local using nzb_tmp as a temporary global index array
822
823!
824!--    nzb_s_outer:
825!--    extend nzb_local east-/westwards first, then north-/southwards
826       nzb_tmp = nzb_local(-nbgp:ny+nbgp,-nbgp:nx+nbgp)
827       DO  j = -1, ny + 1
828          DO  i = 0, nx
829             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i),             &
830                                 nzb_local(j,i+1) )
831          ENDDO
832       ENDDO
833       DO  i = nxl, nxr
834          DO  j = nys, nyn
835             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
836                                     nzb_tmp(j+1,i) )
837          ENDDO
838!
839!--       non-cyclic boundary conditions (overwritten by call of
840!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
841          IF ( nys == 0 )  THEN
842             j = -1
843             nzb_s_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
844          ENDIF
845          IF ( nyn == ny )  THEN
846             j = ny + 1
847             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
848          ENDIF
849       ENDDO
850!
851!--    nzb_w_outer:
852!--    identical to nzb_s_outer
853       nzb_w_outer = nzb_s_outer
854
855!
856!--    nzb_u_inner:
857!--    extend nzb_local rightwards only
858       nzb_tmp = nzb_local(-nbgp:ny+nbgp,-nbgp:nx+nbgp)
859       DO  j = -1, ny + 1
860          DO  i = 0, nx + 1
861             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i) )
862          ENDDO
863       ENDDO
864       nzb_u_inner = nzb_tmp(nysg:nyng,nxlg:nxrg)
865
866!
867!--    nzb_u_outer:
868!--    extend current nzb_tmp (nzb_u_inner) north-/southwards
869       DO  i = nxl, nxr
870          DO  j = nys, nyn
871             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
872                                     nzb_tmp(j+1,i) )
873          ENDDO
874!
875!--       non-cyclic boundary conditions (overwritten by call of
876!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
877          IF ( nys == 0 )  THEN
878             j = -1
879             nzb_u_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
880          ENDIF
881          IF ( nyn == ny )  THEN
882             j = ny + 1
883             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
884          ENDIF
885       ENDDO
886
887!
888!--    nzb_v_inner:
889!--    extend nzb_local northwards only
890       nzb_tmp = nzb_local(-nbgp:ny+nbgp,-nbgp:nx+nbgp)
891       DO  i = -1, nx + 1
892          DO  j = 0, ny + 1
893             nzb_tmp(j,i) = MAX( nzb_local(j-1,i), nzb_local(j,i) )
894          ENDDO
895       ENDDO
896       nzb_v_inner = nzb_tmp(nys-nbgp:nyn+nbgp,nxl-nbgp:nxr+nbgp)
897
898!
899!--    nzb_v_outer:
900!--    extend current nzb_tmp (nzb_v_inner) right-/leftwards
901       DO  j = nys, nyn
902          DO  i = nxl, nxr
903             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i),             &
904                                     nzb_tmp(j,i+1) )
905          ENDDO
906!
907!--       non-cyclic boundary conditions (overwritten by call of
908!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
909          IF ( nxl == 0 )  THEN
910             i = -1
911             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i+1), nzb_tmp(j,i) )
912          ENDIF
913          IF ( nxr == nx )  THEN
914             i = nx + 1
915             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i) )
916          ENDIF
917       ENDDO
918#if ! defined ( __check )
919!
920!--    Exchange of lateral boundary values (parallel computers) and cyclic
921!--    boundary conditions, if applicable.
922!--    Since nzb_s_inner and nzb_w_inner are derived directly from nzb_local
923!--    they do not require exchange and are not included here.
924       CALL exchange_horiz_2d_int( nzb_u_inner )
925       CALL exchange_horiz_2d_int( nzb_u_outer )
926       CALL exchange_horiz_2d_int( nzb_v_inner )
927       CALL exchange_horiz_2d_int( nzb_v_outer )
928       CALL exchange_horiz_2d_int( nzb_w_outer )
929       CALL exchange_horiz_2d_int( nzb_s_outer )
930
931!
932!--    Allocate and set the arrays containing the topography height
933       IF ( myid == 0 )  THEN
934
935          ALLOCATE( zu_s_inner(0:nx+1,0:ny+1), zw_w_inner(0:nx+1,0:ny+1) )
936
937          DO  i = 0, nx + 1
938             DO  j = 0, ny + 1
939                zu_s_inner(i,j) = zu(nzb_local(j,i))
940                zw_w_inner(i,j) = zw(nzb_local(j,i))
941             ENDDO
942          ENDDO
943         
944       ENDIF
945!
946!--    Set flag arrays to be used for masking of grid points
947       DO  i = nxlg, nxrg
948          DO  j = nysg, nyng
949             DO  k = nzb, nzt+1
950                IF ( k <= nzb_s_inner(j,i) )  rflags_s_inner(k,j,i) = 0.0_wp
951                IF ( k <= nzb_s_inner(j,i) )  rflags_invers(j,i,k)  = 0.0_wp
952             ENDDO
953          ENDDO
954       ENDDO
955#endif
956    ENDIF
957
958#if ! defined ( __check )
959!
960!-- Preliminary: to be removed after completion of the topography code!
961!-- Set the former default k index arrays nzb_2d
962    nzb_2d      = nzb
963
964!
965!-- Set the individual index arrays which define the k index from which on
966!-- the usual finite difference form (which does not use surface fluxes) is
967!-- applied
968    IF ( constant_flux_layer  .OR.  use_surface_fluxes )  THEN
969       nzb_diff_u         = nzb_u_inner + 2
970       nzb_diff_v         = nzb_v_inner + 2
971       nzb_diff_s_inner   = nzb_s_inner + 2
972       nzb_diff_s_outer   = nzb_s_outer + 2
973    ELSE
974       nzb_diff_u         = nzb_u_inner + 1
975       nzb_diff_v         = nzb_v_inner + 1
976       nzb_diff_s_inner   = nzb_s_inner + 1
977       nzb_diff_s_outer   = nzb_s_outer + 1
978    ENDIF
979
980!
981!-- Calculation of wall switches and factors required by diffusion_u/v.f90 and
982!-- for limitation of near-wall mixing length l_wall further below
983    corner_nl = 0
984    corner_nr = 0
985    corner_sl = 0
986    corner_sr = 0
987    wall_l    = 0
988    wall_n    = 0
989    wall_r    = 0
990    wall_s    = 0
991
992    DO  i = nxl, nxr
993       DO  j = nys, nyn
994!
995!--       u-component
996          IF ( nzb_u_outer(j,i) > nzb_u_outer(j+1,i) )  THEN
997             wall_u(j,i) = 1.0_wp   ! north wall (location of adjacent fluid)
998             fym(j,i)    = 0.0_wp
999             fyp(j,i)    = 1.0_wp
1000          ELSEIF ( nzb_u_outer(j,i) > nzb_u_outer(j-1,i) )  THEN
1001             wall_u(j,i) = 1.0_wp   ! south wall (location of adjacent fluid)
1002             fym(j,i)    = 1.0_wp
1003             fyp(j,i)    = 0.0_wp
1004          ENDIF
1005!
1006!--       v-component
1007          IF ( nzb_v_outer(j,i) > nzb_v_outer(j,i+1) )  THEN
1008             wall_v(j,i) = 1.0_wp   ! rigth wall (location of adjacent fluid)
1009             fxm(j,i)    = 0.0_wp
1010             fxp(j,i)    = 1.0_wp
1011          ELSEIF ( nzb_v_outer(j,i) > nzb_v_outer(j,i-1) )  THEN
1012             wall_v(j,i) = 1.0_wp   ! left wall (location of adjacent fluid)
1013             fxm(j,i)    = 1.0_wp
1014             fxp(j,i)    = 0.0_wp
1015          ENDIF
1016!
1017!--       w-component, also used for scalars, separate arrays for shear
1018!--       production of tke
1019          IF ( nzb_w_outer(j,i) > nzb_w_outer(j+1,i) )  THEN
1020             wall_e_y(j,i) =  1.0_wp   ! north wall (location of adjacent fluid)
1021             wall_w_y(j,i) =  1.0_wp
1022             fwym(j,i)     =  0.0_wp
1023             fwyp(j,i)     =  1.0_wp
1024          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j-1,i) )  THEN
1025             wall_e_y(j,i) = -1.0_wp   ! south wall (location of adjacent fluid)
1026             wall_w_y(j,i) =  1.0_wp
1027             fwym(j,i)     =  1.0_wp
1028             fwyp(j,i)     =  0.0_wp
1029          ENDIF
1030          IF ( nzb_w_outer(j,i) > nzb_w_outer(j,i+1) )  THEN
1031             wall_e_x(j,i) =  1.0_wp   ! right wall (location of adjacent fluid)
1032             wall_w_x(j,i) =  1.0_wp
1033             fwxm(j,i)     =  0.0_wp
1034             fwxp(j,i)     =  1.0_wp
1035          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j,i-1) )  THEN
1036             wall_e_x(j,i) = -1.0_wp   ! left wall (location of adjacent fluid)
1037             wall_w_x(j,i) =  1.0_wp
1038             fwxm(j,i)     =  1.0_wp
1039             fwxp(j,i)     =  0.0_wp
1040          ENDIF
1041!
1042!--       Wall and corner locations inside buildings for limitation of
1043!--       near-wall mixing length l_wall
1044          IF ( nzb_s_inner(j,i) > nzb_s_inner(j+1,i) )  THEN
1045
1046             wall_n(j,i) = nzb_s_inner(j+1,i) + 1            ! North wall
1047
1048             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1049                corner_nl(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northleft corner
1050                                      nzb_s_inner(j,i-1) ) + 1
1051             ENDIF
1052
1053             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1054                corner_nr(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northright corner
1055                                      nzb_s_inner(j,i+1) ) + 1
1056             ENDIF
1057
1058          ENDIF
1059
1060          IF ( nzb_s_inner(j,i) > nzb_s_inner(j-1,i) )  THEN
1061
1062             wall_s(j,i) = nzb_s_inner(j-1,i) + 1            ! South wall
1063             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1064                corner_sl(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southleft corner
1065                                      nzb_s_inner(j,i-1) ) + 1
1066             ENDIF
1067
1068             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1069                corner_sr(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southright corner
1070                                      nzb_s_inner(j,i+1) ) + 1
1071             ENDIF
1072
1073          ENDIF
1074
1075          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1076             wall_l(j,i) = nzb_s_inner(j,i-1) + 1            ! Left wall
1077          ENDIF
1078
1079          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1080             wall_r(j,i) = nzb_s_inner(j,i+1) + 1            ! Right wall
1081          ENDIF
1082
1083       ENDDO
1084    ENDDO
1085
1086!
1087!-- Calculate wall flag arrays for the multigrid method
1088    IF ( psolver(1:9) == 'multigrid' )  THEN
1089!
1090!--    Gridpoint increment of the current level
1091       inc = 1
1092
1093       DO  l = maximum_grid_level, 1 , -1
1094
1095          nxl_l = nxl_mg(l)
1096          nxr_l = nxr_mg(l)
1097          nys_l = nys_mg(l)
1098          nyn_l = nyn_mg(l)
1099          nzt_l = nzt_mg(l)
1100
1101!
1102!--       Assign the flag level to be calculated
1103          SELECT CASE ( l )
1104             CASE ( 1 )
1105                flags => wall_flags_1
1106             CASE ( 2 )
1107                flags => wall_flags_2
1108             CASE ( 3 )
1109                flags => wall_flags_3
1110             CASE ( 4 )
1111                flags => wall_flags_4
1112             CASE ( 5 )
1113                flags => wall_flags_5
1114             CASE ( 6 )
1115                flags => wall_flags_6
1116             CASE ( 7 )
1117                flags => wall_flags_7
1118             CASE ( 8 )
1119                flags => wall_flags_8
1120             CASE ( 9 )
1121                flags => wall_flags_9
1122             CASE ( 10 )
1123                flags => wall_flags_10
1124          END SELECT
1125
1126!
1127!--       Depending on the grid level, set the respective bits in case of
1128!--       neighbouring walls
1129!--       Bit 0:  wall to the bottom
1130!--       Bit 1:  wall to the top (not realized in remaining PALM code so far)
1131!--       Bit 2:  wall to the south
1132!--       Bit 3:  wall to the north
1133!--       Bit 4:  wall to the left
1134!--       Bit 5:  wall to the right
1135!--       Bit 6:  inside building
1136
1137          flags = 0
1138
1139!
1140!--       In case of masking method, flags are not set and multigrid method
1141!--       works like FFT-solver
1142          IF ( .NOT. masking_method )  THEN
1143
1144             DO  i = nxl_l-1, nxr_l+1
1145                DO  j = nys_l-1, nyn_l+1
1146                   DO  k = nzb, nzt_l+1
1147                         
1148!
1149!--                   Inside/outside building (inside building does not need
1150!--                   further tests for walls)
1151                      IF ( k*inc <= nzb_local(j*inc,i*inc) )  THEN
1152
1153                         flags(k,j,i) = IBSET( flags(k,j,i), 6 )
1154
1155                      ELSE
1156!
1157!--                      Bottom wall
1158                         IF ( (k-1)*inc <= nzb_local(j*inc,i*inc) )  THEN
1159                            flags(k,j,i) = IBSET( flags(k,j,i), 0 )
1160                         ENDIF
1161!
1162!--                      South wall
1163                         IF ( k*inc <= nzb_local((j-1)*inc,i*inc) )  THEN
1164                            flags(k,j,i) = IBSET( flags(k,j,i), 2 )
1165                         ENDIF
1166!
1167!--                      North wall
1168                         IF ( k*inc <= nzb_local((j+1)*inc,i*inc) )  THEN
1169                            flags(k,j,i) = IBSET( flags(k,j,i), 3 )
1170                         ENDIF
1171!
1172!--                      Left wall
1173                         IF ( k*inc <= nzb_local(j*inc,(i-1)*inc) )  THEN
1174                            flags(k,j,i) = IBSET( flags(k,j,i), 4 )
1175                         ENDIF
1176!
1177!--                      Right wall
1178                         IF ( k*inc <= nzb_local(j*inc,(i+1)*inc) )  THEN
1179                            flags(k,j,i) = IBSET( flags(k,j,i), 5 )
1180                         ENDIF
1181
1182                      ENDIF
1183                           
1184                   ENDDO
1185                ENDDO
1186             ENDDO
1187
1188          ENDIF
1189
1190!
1191!--       Test output of flag arrays
1192!          i = nxl_l
1193!          WRITE (9,*)  ' '
1194!          WRITE (9,*)  '*** mg level ', l, ' ***', mg_switch_to_pe0_level
1195!          WRITE (9,*)  '    inc=', inc, '  i =', nxl_l
1196!          WRITE (9,*)  '    nxl_l',nxl_l,' nxr_l=',nxr_l,' nys_l=',nys_l,' nyn_l=',nyn_l
1197!          DO  k = nzt_l+1, nzb, -1
1198!             WRITE (9,'(194(1X,I2))')  ( flags(k,j,i), j = nys_l-1, nyn_l+1 )
1199!          ENDDO
1200
1201          inc = inc * 2
1202
1203       ENDDO
1204
1205    ENDIF
1206!
1207!-- Allocate flags needed for masking walls.
1208    ALLOCATE( wall_flags_0(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
1209              wall_flags_00(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
1210    wall_flags_0  = 0
1211    wall_flags_00 = 0
1212
1213    IF ( scalar_advec == 'ws-scheme' .OR.                                     &
1214         scalar_advec == 'ws-scheme-mono' )  THEN
1215!
1216!--    Set flags to steer the degradation of the advection scheme in advec_ws
1217!--    near topography, inflow- and outflow boundaries as well as bottom and
1218!--    top of model domain. wall_flags_0 remains zero for all non-prognostic
1219!--    grid points.
1220       DO  i = nxl, nxr
1221          DO  j = nys, nyn
1222             DO  k = nzb_s_inner(j,i)+1, nzt
1223!
1224!--             scalar - x-direction
1225!--             WS1 (0), WS3 (1), WS5 (2)
1226                IF ( k <= nzb_s_inner(j,i+1) .OR.                              &
1227                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1228                           .AND. i == nxl   )    .OR.                          &
1229                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1230                           .AND. i == nxr   ) )                                &
1231                THEN
1232                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 0 )
1233                ELSEIF ( k <= nzb_s_inner(j,i+2) .OR. k <= nzb_s_inner(j,i-1)  &
1234                                                 .OR.                          &
1235                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1236                           .AND. i == nxr-1 )    .OR.                          &
1237                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1238                           .AND. i == nxlu  ) )                                &
1239                THEN
1240                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 1 )
1241                ELSE
1242                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 2 )
1243                ENDIF
1244!
1245!--             scalar - y-direction
1246!--             WS1 (3), WS3 (4), WS5 (5)
1247                IF ( k <= nzb_s_inner(j+1,i) .OR.                              &
1248                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1249                           .AND. j == nys   )    .OR.                          &
1250                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1251                           .AND. j == nyn   ) )                                &
1252                THEN
1253                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 3 )
1254!
1255!--             WS3
1256                ELSEIF ( k <= nzb_s_inner(j+2,i) .OR. k <= nzb_s_inner(j-1,i)  &
1257                                                 .OR.                          &
1258                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1259                           .AND. j == nysv  )    .OR.                          &
1260                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1261                           .AND. j == nyn-1 ) )                                &
1262                THEN
1263                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 4 )
1264!
1265!--             WS5
1266                ELSE
1267                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 5 )
1268                ENDIF
1269!
1270!--             scalar - z-direction
1271!--             WS1 (6), WS3 (7), WS5 (8)
1272                flag_set = .FALSE.
1273                IF ( k == nzb_s_inner(j,i) + 1 .OR. k == nzt )  THEN
1274                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 6 )
1275                   flag_set = .TRUE.
1276                ELSEIF ( k == nzb_s_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1277                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 7 )
1278                   flag_set = .TRUE.
1279                ELSEIF ( k > nzb_s_inner(j,i) .AND. .NOT. flag_set )  THEN
1280                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 8 )
1281                ENDIF
1282             ENDDO
1283          ENDDO
1284       ENDDO
1285    ENDIF
1286
1287    IF ( momentum_advec == 'ws-scheme' )  THEN
1288!
1289!--    Set wall_flags_0 to steer the degradation of the advection scheme in advec_ws
1290!--    near topography, inflow- and outflow boundaries as well as bottom and
1291!--    top of model domain. wall_flags_0 remains zero for all non-prognostic
1292!--    grid points.
1293       DO  i = nxl, nxr
1294          DO  j = nys, nyn
1295             DO  k = nzb+1, nzt
1296!
1297!--             At first, set flags to WS1.
1298!--             Since fluxes are swapped in advec_ws.f90, this is necessary to
1299!--             in order to handle the left/south flux.
1300!--             near vertical walls.
1301                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 9 )
1302                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 12 )
1303!
1304!--             u component - x-direction
1305!--             WS1 (9), WS3 (10), WS5 (11)
1306                IF ( k <= nzb_u_inner(j,i+1) .OR.                              &
1307                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1308                           .AND. i <= nxlu  )    .OR.                          &
1309                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1310                           .AND. i == nxr   ) )                                &
1311                THEN
1312                    wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 9 )
1313                ELSEIF ( k <= nzb_u_inner(j,i+2) .OR. k <= nzb_u_inner(j,i-1)  &
1314                                                 .OR.                          &
1315                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1316                           .AND. i == nxr-1 )    .OR.                          &
1317                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1318                           .AND. i == nxlu+1) )                                &
1319                THEN
1320                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 10 )
1321!
1322!--                Clear flag for WS1
1323                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 9 )
1324                ELSE
1325                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 11 )
1326!
1327!--                Clear flag for WS1
1328                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 9 )
1329                ENDIF
1330!
1331!--             u component - y-direction
1332!--             WS1 (12), WS3 (13), WS5 (14)
1333                IF ( k <= nzb_u_inner(j+1,i) .OR.                              &
1334                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1335                           .AND. j == nys   )    .OR.                          &
1336                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1337                           .AND. j == nyn   ) )                                &
1338                THEN
1339                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 12 )
1340                ELSEIF ( k <= nzb_u_inner(j+2,i) .OR. k <= nzb_u_inner(j-1,i)  &
1341                                                 .OR.                          &
1342                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1343                           .AND. j == nysv  )    .OR.                          &
1344                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1345                           .AND. j == nyn-1 ) )                                &
1346                THEN
1347                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 13 )
1348!
1349!--                Clear flag for WS1
1350                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 12 )
1351                ELSE
1352                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 14 )
1353!
1354!--                Clear flag for WS1
1355                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 12 )
1356                ENDIF
1357!
1358!--             u component - z-direction
1359!--             WS1 (15), WS3 (16), WS5 (17)
1360                flag_set = .FALSE.
1361                IF ( k == nzb_u_inner(j,i) + 1 .OR. k == nzt )  THEN
1362                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 15 )
1363                   flag_set = .TRUE.
1364                ELSEIF ( k == nzb_u_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1365                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 16 )
1366                   flag_set = .TRUE.
1367                ELSEIF ( k > nzb_u_inner(j,i) .AND. .NOT. flag_set )  THEN
1368                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 17 )
1369                ENDIF
1370
1371             ENDDO
1372          ENDDO
1373       ENDDO
1374
1375       DO  i = nxl, nxr
1376          DO  j = nys, nyn
1377             DO  k = nzb+1, nzt
1378!
1379!--             At first, set flags to WS1.
1380!--             Since fluxes are swapped in advec_ws.f90, this is necessary to
1381!--             in order to handle the left/south flux.
1382                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 18 )
1383                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 21 )
1384!
1385!--             v component - x-direction
1386!--             WS1 (18), WS3 (19), WS5 (20)
1387                IF ( k <= nzb_v_inner(j,i+1) .OR.                              &
1388                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1389                           .AND. i == nxl   )    .OR.                          &
1390                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1391                           .AND. i == nxr   ) )                                &
1392                THEN
1393                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 18 )
1394!
1395!--             WS3
1396                ELSEIF ( k <= nzb_v_inner(j,i+2) .OR. k <= nzb_v_inner(j,i-1)  &
1397                                                 .OR.                          &
1398                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1399                           .AND. i == nxr-1 )    .OR.                          &
1400                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1401                           .AND. i == nxlu  ) )                                &
1402                THEN
1403                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 19 )
1404!
1405!--                Clear flag for WS1
1406                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 18 )
1407                ELSE
1408                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 20 )
1409!
1410!--                Clear flag for WS1
1411                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 18 )
1412                ENDIF
1413!
1414!--             v component - y-direction
1415!--             WS1 (21), WS3 (22), WS5 (23)
1416                IF ( k <= nzb_v_inner(j+1,i) .OR.                              &
1417                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1418                           .AND. j <= nysv  )    .OR.                          &
1419                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1420                           .AND. j == nyn   ) )                                &
1421                THEN
1422                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 21 )
1423                ELSEIF ( k <= nzb_v_inner(j+2,i) .OR. k <= nzb_v_inner(j-1,i)  &
1424                                                 .OR.                          &
1425                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1426                           .AND. j == nysv+1)    .OR.                          &
1427                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1428                           .AND. j == nyn-1 ) )                                &
1429                THEN
1430                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 22 )
1431!
1432!--                Clear flag for WS1
1433                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 21 )
1434                ELSE
1435                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 23 )
1436!
1437!--                Clear flag for WS1
1438                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 21 )
1439                ENDIF
1440!
1441!--             v component - z-direction
1442!--             WS1 (24), WS3 (25), WS5 (26)
1443                flag_set = .FALSE.
1444                IF ( k == nzb_v_inner(j,i) + 1 .OR. k == nzt )  THEN
1445                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 24 )
1446                   flag_set = .TRUE.
1447                ELSEIF ( k == nzb_v_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1448                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 25 )
1449                   flag_set = .TRUE.
1450                ELSEIF ( k > nzb_v_inner(j,i) .AND. .NOT. flag_set )  THEN
1451                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 26 )
1452                ENDIF
1453
1454             ENDDO
1455          ENDDO
1456       ENDDO
1457       DO  i = nxl, nxr
1458          DO  j = nys, nyn
1459             DO  k = nzb+1, nzt
1460!
1461!--             At first, set flags to WS1.
1462!--             Since fluxes are swapped in advec_ws.f90, this is necessary to
1463!--             in order to handle the left/south flux.
1464                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 27 )
1465                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 30 )
1466!
1467!--             w component - x-direction
1468!--             WS1 (27), WS3 (28), WS5 (29)
1469                IF ( k <= nzb_w_inner(j,i+1) .OR.                              &
1470                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1471                           .AND. i == nxl   )    .OR.                          &
1472                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1473                           .AND. i == nxr   ) )                                &
1474                THEN
1475                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 27 )
1476                ELSEIF ( k <= nzb_w_inner(j,i+2) .OR. k <= nzb_w_inner(j,i-1)  &
1477                                                 .OR.                          &
1478                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1479                           .AND. i == nxr-1 )    .OR.                          &
1480                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1481                           .AND. i == nxlu  ) )                                &
1482                THEN
1483                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 28 )
1484!
1485!--                Clear flag for WS1
1486                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 27 )
1487                ELSE
1488                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i),29 )
1489!
1490!--                Clear flag for WS1
1491                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 27 )
1492                ENDIF
1493!
1494!--             w component - y-direction
1495!--             WS1 (30), WS3 (31), WS5 (32)
1496                IF ( k <= nzb_w_inner(j+1,i) .OR.                              &
1497                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1498                           .AND. j == nys   )    .OR.                          &
1499                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1500                           .AND. j == nyn   ) )                                &
1501                THEN
1502                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 30 )
1503                ELSEIF ( k <= nzb_w_inner(j+2,i) .OR. k <= nzb_w_inner(j-1,i)  &
1504                                                 .OR.                          &
1505                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1506                           .AND. j == nysv  )    .OR.                          &
1507                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1508                           .AND. j == nyn-1 ) )                                &
1509                THEN
1510                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 31 )
1511!
1512!--                Clear flag for WS1
1513                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 30 )
1514                ELSE
1515                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 0 )
1516!
1517!--                Clear flag for WS1
1518                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 30 )
1519                ENDIF
1520!
1521!--             w component - z-direction
1522!--             WS1 (33), WS3 (34), WS5 (35)
1523                flag_set = .FALSE.
1524                IF ( k == nzb_w_inner(j,i) .OR. k == nzb_w_inner(j,i) + 1      &
1525                                           .OR. k == nzt )  THEN
1526!
1527!--                Please note, at k == nzb_w_inner(j,i) a flag is explictely
1528!--                set, although this is not a prognostic level. However,
1529!--                contrary to the advection of u,v and s this is necessary
1530!--                because flux_t(nzb_w_inner(j,i)) is used for the tendency
1531!--                at k == nzb_w_inner(j,i)+1.
1532                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 1 )
1533                   flag_set = .TRUE.
1534                ELSEIF ( k == nzb_w_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1535                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 2 )
1536                   flag_set = .TRUE.
1537                ELSEIF ( k > nzb_w_inner(j,i) .AND. .NOT. flag_set )  THEN
1538                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 3 )
1539                ENDIF
1540
1541             ENDDO
1542          ENDDO
1543       ENDDO
1544
1545    ENDIF
1546
1547!
1548!-- Exchange 3D integer wall_flags.
1549    IF ( momentum_advec == 'ws-scheme' .OR. scalar_advec == 'ws-scheme'     &
1550    .OR. scalar_advec == 'ws-scheme-mono' )  THEN 
1551!
1552!--    Exchange ghost points for advection flags
1553       CALL exchange_horiz_int( wall_flags_0,  nbgp )
1554       CALL exchange_horiz_int( wall_flags_00, nbgp )
1555!
1556!--    Set boundary flags at inflow and outflow boundary in case of
1557!--    non-cyclic boundary conditions.
1558       IF ( inflow_l .OR. outflow_l .OR. nest_bound_l )  THEN
1559          wall_flags_0(:,:,nxl-1)  = wall_flags_0(:,:,nxl)
1560          wall_flags_00(:,:,nxl-1) = wall_flags_00(:,:,nxl)
1561       ENDIF
1562
1563       IF ( inflow_r .OR. outflow_r .OR. nest_bound_r )  THEN
1564          wall_flags_0(:,:,nxr+1)  = wall_flags_0(:,:,nxr)
1565          wall_flags_00(:,:,nxr+1) = wall_flags_00(:,:,nxr)
1566       ENDIF
1567
1568       IF ( inflow_n .OR. outflow_n .OR. nest_bound_n )  THEN
1569          wall_flags_0(:,nyn+1,:)  = wall_flags_0(:,nyn,:)
1570          wall_flags_00(:,nyn+1,:) = wall_flags_00(:,nyn,:)
1571       ENDIF
1572
1573       IF ( inflow_s .OR. outflow_s  .OR. nest_bound_s )  THEN
1574          wall_flags_0(:,nys-1,:)  = wall_flags_0(:,nys,:)
1575          wall_flags_00(:,nys-1,:) = wall_flags_00(:,nys,:)
1576       ENDIF
1577
1578    ENDIF
1579
1580!
1581!-- In case of topography: limit near-wall mixing length l_wall further:
1582!-- Go through all points of the subdomain one by one and look for the closest
1583!-- surface
1584    IF ( TRIM(topography) /= 'flat' )  THEN
1585       DO  i = nxl, nxr
1586          DO  j = nys, nyn
1587
1588             nzb_si = nzb_s_inner(j,i)
1589             vi     = vertical_influence(nzb_si)
1590
1591             IF ( wall_n(j,i) > 0 )  THEN
1592!
1593!--             North wall (y distance)
1594                DO  k = wall_n(j,i), nzb_si
1595                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i), 0.5_wp * dy )
1596                ENDDO
1597!
1598!--             Above North wall (yz distance)
1599                DO  k = nzb_si + 1, nzb_si + vi
1600                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i),                     &
1601                                          SQRT( 0.25_wp * dy**2 +              &
1602                                          ( zu(k) - zw(nzb_si) )**2 ) )
1603                ENDDO
1604!
1605!--             Northleft corner (xy distance)
1606                IF ( corner_nl(j,i) > 0 )  THEN
1607                   DO  k = corner_nl(j,i), nzb_si
1608                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1), &
1609                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1610                   ENDDO
1611!
1612!--                Above Northleft corner (xyz distance)
1613                   DO  k = nzb_si + 1, nzb_si + vi
1614                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1),              &
1615                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1616                                            ( zu(k) - zw(nzb_si) )**2 ) )
1617                   ENDDO
1618                ENDIF
1619!
1620!--             Northright corner (xy distance)
1621                IF ( corner_nr(j,i) > 0 )  THEN
1622                   DO  k = corner_nr(j,i), nzb_si
1623                       l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1),             &
1624                                                0.5_wp * SQRT( dx**2 + dy**2 ) )
1625                   ENDDO
1626!
1627!--                Above northright corner (xyz distance)
1628                   DO  k = nzb_si + 1, nzb_si + vi
1629                      l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1),              &
1630                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1631                                            ( zu(k) - zw(nzb_si) )**2 ) )
1632                   ENDDO
1633                ENDIF
1634             ENDIF
1635
1636             IF ( wall_s(j,i) > 0 )  THEN
1637!
1638!--             South wall (y distance)
1639                DO  k = wall_s(j,i), nzb_si
1640                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i), 0.5_wp * dy )
1641                ENDDO
1642!
1643!--             Above south wall (yz distance)
1644                DO  k = nzb_si + 1, nzb_si + vi
1645                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i),                     &
1646                                          SQRT( 0.25_wp * dy**2 +              &
1647                                          ( zu(k) - zw(nzb_si) )**2 ) )
1648                ENDDO
1649!
1650!--             Southleft corner (xy distance)
1651                IF ( corner_sl(j,i) > 0 )  THEN
1652                   DO  k = corner_sl(j,i), nzb_si
1653                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1654                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1655                   ENDDO
1656!
1657!--                Above southleft corner (xyz distance)
1658                   DO  k = nzb_si + 1, nzb_si + vi
1659                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1660                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1661                                            ( zu(k) - zw(nzb_si) )**2 ) )
1662                   ENDDO
1663                ENDIF
1664!
1665!--             Southright corner (xy distance)
1666                IF ( corner_sr(j,i) > 0 )  THEN
1667                   DO  k = corner_sr(j,i), nzb_si
1668                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1669                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1670                   ENDDO
1671!
1672!--                Above southright corner (xyz distance)
1673                   DO  k = nzb_si + 1, nzb_si + vi
1674                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1675                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1676                                            ( zu(k) - zw(nzb_si) )**2 ) )
1677                   ENDDO
1678                ENDIF
1679
1680             ENDIF
1681
1682             IF ( wall_l(j,i) > 0 )  THEN
1683!
1684!--             Left wall (x distance)
1685                DO  k = wall_l(j,i), nzb_si
1686                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1), 0.5_wp * dx )
1687                ENDDO
1688!
1689!--             Above left wall (xz distance)
1690                DO  k = nzb_si + 1, nzb_si + vi
1691                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1),                     &
1692                                       SQRT( 0.25_wp * dx**2 +                 &
1693                                       ( zu(k) - zw(nzb_si) )**2 ) )
1694                ENDDO
1695             ENDIF
1696
1697             IF ( wall_r(j,i) > 0 )  THEN
1698!
1699!--             Right wall (x distance)
1700                DO  k = wall_r(j,i), nzb_si
1701                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1), 0.5_wp * dx )
1702                ENDDO
1703!
1704!--             Above right wall (xz distance)
1705                DO  k = nzb_si + 1, nzb_si + vi
1706                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1),                     &
1707                                          SQRT( 0.25_wp * dx**2 +              &
1708                                          ( zu(k) - zw(nzb_si) )**2 ) )
1709                ENDDO
1710
1711             ENDIF
1712
1713          ENDDO
1714       ENDDO
1715
1716    ENDIF
1717
1718!
1719!-- Multiplication with wall_adjustment_factor
1720    l_wall = wall_adjustment_factor * l_wall
1721
1722!
1723!-- Set lateral boundary conditions for l_wall
1724    CALL exchange_horiz( l_wall, nbgp )
1725
1726    DEALLOCATE( corner_nl, corner_nr, corner_sl, corner_sr, nzb_local, &
1727                nzb_tmp, vertical_influence, wall_l, wall_n, wall_r, wall_s )
1728
1729#endif
1730
1731 END SUBROUTINE init_grid
Note: See TracBrowser for help on using the repository browser.