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

Last change on this file since 1942 was 1942, checked in by suehring, 8 years ago

filter topography by filling holes resolved by only one grid point ;initialization of wall_flags_0 and wall_flags_00 moved to advec_ws

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