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

Last change on this file since 2069 was 2038, checked in by knoop, 8 years ago

last commit documented

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