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

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

Bugfix in definition of generic topography; missing check for microphysics added

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