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

Last change on this file since 1875 was 1846, checked in by raasch, 9 years ago

last commit documented

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