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

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

last commit documented

  • Property svn:keywords set to Id
File size: 61.8 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! -----------------
[1969]21!
[1983]22!
[1969]23! Former revisions:
24! -----------------
25! $Id: init_grid.f90 1983 2016-08-01 11:06:41Z 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
[240]644       CASE ( 'single_street_canyon' )
645!
646!--       Single quasi-2D street canyon of infinite length in x or y direction.
647!--       The canyon is centered in the other direction by default.
[1322]648          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[240]649!
650!--          Street canyon in y direction
651             cwx = NINT( canyon_width_x / dx )
[1322]652             IF ( canyon_wall_left == 9999999.9_wp )  THEN
[240]653                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
654             ENDIF
655             cxl = NINT( canyon_wall_left / dx )
656             cxr = cxl + cwx
657
[1322]658          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[240]659!
660!--          Street canyon in x direction
661             cwy = NINT( canyon_width_y / dy )
[1322]662             IF ( canyon_wall_south == 9999999.9_wp )  THEN
[240]663                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
664             ENDIF
665             cys = NINT( canyon_wall_south / dy )
666             cyn = cys + cwy
667
668          ELSE
[254]669             
670             message_string = 'no street canyon width given'
671             CALL message( 'init_grid', 'PA0204', 1, 2, 0, 6, 0 )
672 
[240]673          ENDIF
674
[1675]675          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
676          IF ( ABS( zw(ch  ) - canyon_height ) == &
677               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
678
[240]679          dp_level_ind_b = ch
680!
681!--       Street canyon size has to meet some requirements
[1322]682          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[1353]683             IF ( ( cxl < 1 ) .OR. ( cxr > nx-1 ) .OR. ( cwx < 3 ) .OR.        &
[240]684               ( ch < 3 ) )  THEN
[1353]685                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
686                                           '&cxl=', cxl, 'cxr=', cxr,          &
687                                           'cwx=', cwx,                        &
[254]688                                           'ch=', ch, 'nx=', nx, 'ny=', ny
689                CALL message( 'init_grid', 'PA0205', 1, 2, 0, 6, 0 ) 
[240]690             ENDIF
[1322]691          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[1353]692             IF ( ( cys < 1 ) .OR. ( cyn > ny-1 ) .OR. ( cwy < 3 ) .OR.        &
[240]693               ( ch < 3 ) )  THEN
[1353]694                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
695                                           '&cys=', cys, 'cyn=', cyn,          &
696                                           'cwy=', cwy,                        &
[254]697                                           'ch=', ch, 'nx=', nx, 'ny=', ny
698                CALL message( 'init_grid', 'PA0206', 1, 2, 0, 6, 0 ) 
[240]699             ENDIF
700          ENDIF
[1353]701          IF ( canyon_width_x /= 9999999.9_wp .AND.                            &                 
702               canyon_width_y /= 9999999.9_wp )  THEN
703             message_string = 'inconsistent canyon parameters:' //             &   
704                              '&street canyon can only be oriented' //         &
[254]705                              '&either in x- or in y-direction'
706             CALL message( 'init_grid', 'PA0207', 1, 2, 0, 6, 0 )
[240]707          ENDIF
708
709          nzb_local = ch
[1322]710          IF ( canyon_width_x /= 9999999.9_wp )  THEN
[1968]711             IF ( cxl <= nxr  .AND.  cxr >= nxl )                              &
712                nzb_local(:,MAX(nxl,cxl+1):MIN(nxr,cxr-1)) = 0
[1322]713          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
[1968]714             IF ( cys <= nyn  .AND.  cyn >= nys )                              &         
715                nzb_local(MAX(nys,cys+1):MIN(nyn,cyn-1),:) = 0
[240]716          ENDIF
717
[1]718       CASE ( 'read_from_file' )
[759]719
[1968]720          ALLOCATE ( topo_height(nys:nyn,nxl:nxr) )
[818]721
[759]722          DO  ii = 0, io_blocks-1
723             IF ( ii == io_group )  THEN
724
[1]725!
[759]726!--             Arbitrary irregular topography data in PALM format (exactly
727!--             matching the grid size and total domain size)
[1779]728                OPEN( 90, FILE='TOPOGRAPHY_DATA'//TRIM( coupling_char ),       &
729                          STATUS='OLD', FORM='FORMATTED', ERR=10 )
[1968]730!
731!--             Read topography PE-wise. Rows are read from nyn to nys, columns
732!--             are read from nxl to nxr. At first, ny-nyn rows need to be skipped.
733                skip_n_rows = 0
734                DO WHILE ( skip_n_rows < ny - nyn )
735                   READ( 90, * ) 
736                   skip_n_rows = skip_n_rows + 1
[759]737                ENDDO
[1968]738!
739!--             Read data from nyn to nys and nxl to nxr. Therefore, skip
740!--             column until nxl-1 is reached
741                DO  j = nyn, nys, -1
742                   READ( 90, *, ERR=11, END=11 )                               &
743                                              ( dum, i = 0, nxl-1 ),           &
744                                              ( topo_height(j,i), i = nxl, nxr )
745                ENDDO
[759]746
747                GOTO 12
748         
[1779]749 10             message_string = 'file TOPOGRAPHY'//TRIM( coupling_char )//    &
750                                 ' does not exist'
[759]751                CALL message( 'init_grid', 'PA0208', 1, 2, 0, 6, 0 )
752
[1779]753 11             message_string = 'errors in file TOPOGRAPHY_DATA'//            &
754                                 TRIM( coupling_char )
[759]755                CALL message( 'init_grid', 'PA0209', 1, 2, 0, 6, 0 )
756
757 12             CLOSE( 90 )
758
759             ENDIF
[1804]760#if defined( __parallel )
[759]761             CALL MPI_BARRIER( comm2d, ierr )
762#endif
[559]763          ENDDO
[759]764
[1]765!
766!--       Calculate the index height of the topography
[1968]767          nzb_local = 0
768          DO  i = nxl, nxr
769             DO  j = nys, nyn
[1675]770                nzb_local(j,i) = MINLOC( ABS( zw - topo_height(j,i) ), 1 ) - 1
771                IF ( ABS( zw(nzb_local(j,i)  ) - topo_height(j,i) ) == &
772                     ABS( zw(nzb_local(j,i)+1) - topo_height(j,i) )    )  &
773                   nzb_local(j,i) = nzb_local(j,i) + 1
[1]774             ENDDO
775          ENDDO
[818]776
777          DEALLOCATE ( topo_height )
[1942]778!
779!--       Filter topography, i.e. fill holes resolved by only one grid point. 
780!--       Such holes are suspected to lead to velocity blow-ups as continuity
781!--       equation on discrete grid cannot be fulfilled in such case.
782!--       For now, check only for holes and fill them to the lowest height level
783!--       of the directly adjoining grid points along x- and y- direction.
784!--       Before checking for holes, set lateral boundary conditions for
785!--       topography. After hole-filling, boundary conditions must be set again!
[1968]786          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
787         
788          IF ( .NOT. bc_ns_cyc )  THEN
789             IF ( nys == 0  )  nzb_local(-1,:)   = nzb_local(0,:)
790             IF ( nyn == ny )  nzb_local(ny+1,:) = nzb_local(ny,:)
[1942]791          ENDIF
[1910]792
[1968]793          IF ( .NOT. bc_lr_cyc )  THEN
794             IF ( nxl == 0  )  nzb_local(:,-1)   = nzb_local(:,0)
795             IF ( nxr == nx )  nzb_local(:,nx+1) = nzb_local(:,nx)         
[1942]796          ENDIF
797
[1968]798          num_hole_l = 0
799          DO i = nxl, nxr
800             DO j = nys, nyn
[1942]801
802                num_wall = 0
803
804                IF ( nzb_local(j-1,i) > nzb_local(j,i) )                       &
805                   num_wall = num_wall + 1
806                IF ( nzb_local(j+1,i) > nzb_local(j,i) )                       &
807                   num_wall = num_wall + 1
808                IF ( nzb_local(j,i-1) > nzb_local(j,i) )                       &
809                   num_wall = num_wall + 1
810                IF ( nzb_local(j,i+1) > nzb_local(j,i) )                       &
811                   num_wall = num_wall + 1
812
813                IF ( num_wall == 4 )  THEN
814                   nzb_local(j,i) = MIN( nzb_local(j-1,i), nzb_local(j+1,i),   &
815                                         nzb_local(j,i-1), nzb_local(j,i+1) )
[1968]816                   num_hole_l     = num_hole_l + 1
[1942]817                ENDIF
818             ENDDO
819          ENDDO
[114]820!
[1968]821!--       Count the total number of holes, required for informative message.
822#if defined( __parallel )
823          CALL MPI_ALLREDUCE( num_hole_l, num_hole, 1, MPI_INTEGER, MPI_SUM,   &
824                              comm2d, ierr )
825#else
826          num_hole = num_hole_l
827#endif   
828!
[1942]829!--       Create an informative message if any hole was removed.
[1968]830          IF ( num_hole > 0 )  THEN
[1942]831             WRITE( message_string, * ) num_hole, 'hole(s) resolved by only '//&
832                                                  'one grid point were filled'
833             CALL message( 'init_grid', 'PA0430', 0, 0, 0, 6, 0 )
834          ENDIF
835!
[1968]836!--       Exchange ghost-points, as well as add cyclic or Neumann boundary
837!--       conditions.
838          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
839         
840          IF ( .NOT. bc_ns_cyc )  THEN
841             IF ( nys == 0  )  nzb_local(-1,:)   = nzb_local(0,:)
842             IF ( nyn == ny )  nzb_local(ny+1,:) = nzb_local(ny,:)
[1910]843          ENDIF
[667]844
[1968]845          IF ( .NOT. bc_lr_cyc )  THEN
846             IF ( nxl == 0  )  nzb_local(:,-1)   = nzb_local(:,0)
847             IF ( nxr == nx )  nzb_local(:,nx+1) = nzb_local(:,nx)         
[1910]848          ENDIF
849
[1]850       CASE DEFAULT
851!
852!--       The DEFAULT case is reached either if the parameter topography
[217]853!--       contains a wrong character string or if the user has defined a special
[1]854!--       case in the user interface. There, the subroutine user_init_grid
855!--       checks which of these two conditions applies.
[1968]856          CALL user_init_grid( nzb_local )
[1]857
858    END SELECT
859!
[861]860!-- Determine the maximum level of topography. Furthermore it is used for
861!-- steering the degradation of order of the applied advection scheme.
[978]862!-- In case of non-cyclic lateral boundaries, the order of the advection
[1968]863!-- scheme has to be reduced up to nzt (required at the lateral boundaries).
864#if defined( __parallel )
865    CALL MPI_ALLREDUCE( MAXVAL( nzb_local ) + 1, nzb_max, 1, MPI_INTEGER,      &
866                        MPI_MAX, comm2d, ierr )
867#else
[1677]868    nzb_max = MAXVAL( nzb_local ) + 1
[1968]869#endif
[1353]870    IF ( inflow_l .OR. outflow_l .OR. inflow_r .OR. outflow_r .OR.             &
[1762]871         inflow_n .OR. outflow_n .OR. inflow_s .OR. outflow_s .OR.             &
872         nest_domain )                                                         &
873    THEN
874       nzb_max = nzt
[978]875    ENDIF
876
[861]877!
[1]878!-- Consistency checks and index array initialization are only required for
[217]879!-- non-flat topography, also the initialization of topography height arrays
[49]880!-- zu_s_inner and zw_w_inner
[1]881    IF ( TRIM( topography ) /= 'flat' )  THEN
[1968]882#if defined( __parallel )
883       CALL MPI_ALLREDUCE( MAXVAL( nzb_local ), nzb_local_max, 1, MPI_INTEGER, &
884                           MPI_MAX, comm2d, ierr )
[1982]885       CALL MPI_ALLREDUCE( MINVAL( nzb_local ), nzb_local_min, 1, MPI_INTEGER, &
[1968]886                           MPI_MIN, comm2d, ierr )                           
887#else
888       nzb_local_max = MAXVAL( nzb_local )
889       nzb_local_min = MINVAL( nzb_local )
890#endif
[1982]891
[1]892!
893!--    Consistency checks
[1968]894       IF ( nzb_local_min < 0  .OR.  nzb_local_max  > nz + 1 )  THEN
[1353]895          WRITE( message_string, * ) 'nzb_local values are outside the',       &
896                                'model domain',                                &
[1968]897                                '&MINVAL( nzb_local ) = ', nzb_local_min,      &
898                                '&MAXVAL( nzb_local ) = ', nzb_local_max
[254]899          CALL message( 'init_grid', 'PA0210', 1, 2, 0, 6, 0 )
[1]900       ENDIF
901
[217]902       IF ( topography_grid_convention == 'cell_edge' )  THEN
[134]903!
[217]904!--       The array nzb_local as defined using the 'cell_edge' convention
905!--       describes the actual total size of topography which is defined at the
906!--       cell edges where u=0 on the topography walls in x-direction and v=0
907!--       on the topography walls in y-direction. However, PALM uses individual
908!--       arrays nzb_u|v|w|s_inner|outer that are based on nzb_s_inner.
909!--       Therefore, the extent of topography in nzb_local is now reduced by
910!--       1dx at the E topography walls and by 1dy at the N topography walls
[1968]911!--       to form the basis for nzb_s_inner.
912!--       Note, the reverse memory access (i-j instead of j-i) is absolutely
913!--       required at this point.
914          DO  j = nys+1, nyn+1
915             DO  i = nxl-1, nxr
[217]916                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j,i+1) )
917             ENDDO
[134]918          ENDDO
[1968]919!
920!--       Exchange ghost points
921          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
922
923          DO  i = nxl, nxr+1
924             DO  j = nys-1, nyn
[217]925                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j+1,i) )
926             ENDDO
[134]927          ENDDO
[1968]928!
929!--       Exchange ghost points         
930          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
[217]931       ENDIF
[1]932!
933!--    Initialize index arrays nzb_s_inner and nzb_w_inner
[1968]934       nzb_s_inner = nzb_local
935       nzb_w_inner = nzb_local
[1]936
937!
938!--    Initialize remaining index arrays:
939!--    first pre-initialize them with nzb_s_inner...
940       nzb_u_inner = nzb_s_inner
941       nzb_u_outer = nzb_s_inner
942       nzb_v_inner = nzb_s_inner
943       nzb_v_outer = nzb_s_inner
944       nzb_w_outer = nzb_s_inner
945       nzb_s_outer = nzb_s_inner
946
947!
948!--    ...then extend pre-initialized arrays in their according directions
949!--    based on nzb_local using nzb_tmp as a temporary global index array
950
951!
952!--    nzb_s_outer:
953!--    extend nzb_local east-/westwards first, then north-/southwards
[1968]954       nzb_tmp = nzb_local
955       DO  j = nys, nyn
956          DO  i = nxl, nxr
[1353]957             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i),             &
[1]958                                 nzb_local(j,i+1) )
959          ENDDO
960       ENDDO
[1968]961       
962       CALL exchange_horiz_2d_int( nzb_tmp, nys, nyn, nxl, nxr, nbgp )
963       
[1]964       DO  i = nxl, nxr
965          DO  j = nys, nyn
[1353]966             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
[1]967                                     nzb_tmp(j+1,i) )
968          ENDDO
969!
970!--       non-cyclic boundary conditions (overwritten by call of
971!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
972          IF ( nys == 0 )  THEN
973             j = -1
974             nzb_s_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
975          ENDIF
[1743]976          IF ( nyn == ny )  THEN
[1]977             j = ny + 1
978             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
979          ENDIF
980       ENDDO
981!
982!--    nzb_w_outer:
983!--    identical to nzb_s_outer
984       nzb_w_outer = nzb_s_outer
985
986!
987!--    nzb_u_inner:
988!--    extend nzb_local rightwards only
[1968]989       nzb_tmp = nzb_local
990       DO  j = nys, nyn
991          DO  i = nxl, nxr
[1]992             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i) )
993          ENDDO
994       ENDDO
[1968]995       
996       CALL exchange_horiz_2d_int( nzb_tmp, nys, nyn, nxl, nxr, nbgp )
997       
998       nzb_u_inner = nzb_tmp
[1]999!
1000!--    nzb_u_outer:
1001!--    extend current nzb_tmp (nzb_u_inner) north-/southwards
1002       DO  i = nxl, nxr
1003          DO  j = nys, nyn
[1353]1004             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
[1]1005                                     nzb_tmp(j+1,i) )
1006          ENDDO
1007!
1008!--       non-cyclic boundary conditions (overwritten by call of
1009!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
1010          IF ( nys == 0 )  THEN
1011             j = -1
1012             nzb_u_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
1013          ENDIF
[1743]1014          IF ( nyn == ny )  THEN
[1]1015             j = ny + 1
1016             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
1017          ENDIF
1018       ENDDO
1019
1020!
1021!--    nzb_v_inner:
1022!--    extend nzb_local northwards only
[1968]1023       nzb_tmp = nzb_local
1024       DO  i = nxl, nxr
1025          DO  j = nys, nyn
[1]1026             nzb_tmp(j,i) = MAX( nzb_local(j-1,i), nzb_local(j,i) )
1027          ENDDO
1028       ENDDO
[1968]1029       
1030       CALL exchange_horiz_2d_int( nzb_tmp, nys, nyn, nxl, nxr, nbgp )     
1031       nzb_v_inner = nzb_tmp
[1]1032
1033!
1034!--    nzb_v_outer:
1035!--    extend current nzb_tmp (nzb_v_inner) right-/leftwards
1036       DO  j = nys, nyn
1037          DO  i = nxl, nxr
[1353]1038             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i),             &
[1]1039                                     nzb_tmp(j,i+1) )
1040          ENDDO
1041!
1042!--       non-cyclic boundary conditions (overwritten by call of
1043!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
1044          IF ( nxl == 0 )  THEN
1045             i = -1
1046             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i+1), nzb_tmp(j,i) )
1047          ENDIF
1048          IF ( nxr == nx )  THEN
1049             i = nx + 1
1050             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i) )
1051          ENDIF
1052       ENDDO
[1804]1053
[1]1054!
1055!--    Exchange of lateral boundary values (parallel computers) and cyclic
1056!--    boundary conditions, if applicable.
1057!--    Since nzb_s_inner and nzb_w_inner are derived directly from nzb_local
1058!--    they do not require exchange and are not included here.
[1968]1059       CALL exchange_horiz_2d_int( nzb_u_inner, nys, nyn, nxl, nxr, nbgp )
1060       CALL exchange_horiz_2d_int( nzb_u_outer, nys, nyn, nxl, nxr, nbgp )
1061       CALL exchange_horiz_2d_int( nzb_v_inner, nys, nyn, nxl, nxr, nbgp )
1062       CALL exchange_horiz_2d_int( nzb_v_outer, nys, nyn, nxl, nxr, nbgp )
1063       CALL exchange_horiz_2d_int( nzb_w_outer, nys, nyn, nxl, nxr, nbgp )
1064       CALL exchange_horiz_2d_int( nzb_s_outer, nys, nyn, nxl, nxr, nbgp )
[1]1065
[49]1066!
1067!--    Allocate and set the arrays containing the topography height
[1968]1068       ALLOCATE( zu_s_inner(0:nx+1,0:ny+1), zw_w_inner(0:nx+1,0:ny+1),         &
1069                 zu_s_inner_l(0:nx+1,0:ny+1), zw_w_inner_l(0:nx+1,0:ny+1) )
1070                 
1071       zu_s_inner   = 0.0_wp
1072       zw_w_inner   = 0.0_wp
1073       zu_s_inner_l = 0.0_wp
1074       zw_w_inner_l = 0.0_wp
1075       
1076       DO  i = nxl, nxr
1077          DO  j = nys, nyn
1078             zu_s_inner_l(i,j) = zu(nzb_local(j,i))
1079             zw_w_inner_l(i,j) = zw(nzb_local(j,i))
1080          ENDDO
1081       ENDDO
1082       
1083#if defined( __parallel )
1084       CALL MPI_REDUCE( zu_s_inner_l, zu_s_inner, (nx+2)*(ny+2),         &
1085                           MPI_REAL, MPI_SUM, 0, comm2d, ierr )       
1086       CALL MPI_REDUCE( zw_w_inner_l, zw_w_inner, (nx+2)*(ny+2),         &
1087                           MPI_REAL, MPI_SUM, 0, comm2d, ierr ) 
1088#else
1089       zu_s_inner = zu_s_inner_l
1090       zw_w_inner = zw_w_inner_l
1091#endif
[49]1092
[1968]1093      DEALLOCATE( zu_s_inner_l, zw_w_inner_l )
1094      IF ( myid /= 0 )  DEALLOCATE( zu_s_inner, zw_w_inner )
[1221]1095!
[1968]1096!--   Set south and left ghost points, required for netcdf output
1097      IF ( myid == 0 )  THEN
1098         IF( bc_lr_cyc )  THEN
1099            zu_s_inner(nx+1,:) = zu_s_inner(0,:)
1100            zw_w_inner(nx+1,:) = zw_w_inner(0,:)
1101         ELSE
1102            zu_s_inner(nx+1,:) = zu_s_inner(nx,:)
1103            zw_w_inner(nx+1,:) = zw_w_inner(nx,:)
1104         ENDIF
1105         IF( bc_ns_cyc )  THEN
1106            zu_s_inner(:,ny+1) = zu_s_inner(:,0)
1107            zw_w_inner(:,ny+1) = zw_w_inner(:,0)
1108         ELSE
1109            zu_s_inner(:,ny+1) = zu_s_inner(:,ny)
1110            zw_w_inner(:,ny+1) = zw_w_inner(:,ny)
1111         ENDIF
1112      ENDIF
1113!
[1221]1114!--    Set flag arrays to be used for masking of grid points
1115       DO  i = nxlg, nxrg
1116          DO  j = nysg, nyng
1117             DO  k = nzb, nzt+1
[1353]1118                IF ( k <= nzb_s_inner(j,i) )  rflags_s_inner(k,j,i) = 0.0_wp
1119                IF ( k <= nzb_s_inner(j,i) )  rflags_invers(j,i,k)  = 0.0_wp
[1221]1120             ENDDO
1121          ENDDO
1122       ENDDO
[1804]1123
[1]1124    ENDIF
[1968]1125!
1126!-- Deallocate temporary array, as it might be reused for different
1127!-- grid-levels further below.
1128    DEALLOCATE( nzb_tmp )
[1]1129
1130!
1131!-- Set the individual index arrays which define the k index from which on
1132!-- the usual finite difference form (which does not use surface fluxes) is
1133!-- applied
[1691]1134    IF ( constant_flux_layer  .OR.  use_surface_fluxes )  THEN
[1]1135       nzb_diff_u         = nzb_u_inner + 2
1136       nzb_diff_v         = nzb_v_inner + 2
1137       nzb_diff_s_inner   = nzb_s_inner + 2
1138       nzb_diff_s_outer   = nzb_s_outer + 2
1139    ELSE
1140       nzb_diff_u         = nzb_u_inner + 1
1141       nzb_diff_v         = nzb_v_inner + 1
1142       nzb_diff_s_inner   = nzb_s_inner + 1
1143       nzb_diff_s_outer   = nzb_s_outer + 1
1144    ENDIF
1145
1146!
1147!-- Calculation of wall switches and factors required by diffusion_u/v.f90 and
1148!-- for limitation of near-wall mixing length l_wall further below
1149    corner_nl = 0
1150    corner_nr = 0
1151    corner_sl = 0
1152    corner_sr = 0
1153    wall_l    = 0
1154    wall_n    = 0
1155    wall_r    = 0
1156    wall_s    = 0
1157
1158    DO  i = nxl, nxr
1159       DO  j = nys, nyn
1160!
1161!--       u-component
1162          IF ( nzb_u_outer(j,i) > nzb_u_outer(j+1,i) )  THEN
[1353]1163             wall_u(j,i) = 1.0_wp   ! north wall (location of adjacent fluid)
1164             fym(j,i)    = 0.0_wp
1165             fyp(j,i)    = 1.0_wp
[1]1166          ELSEIF ( nzb_u_outer(j,i) > nzb_u_outer(j-1,i) )  THEN
[1353]1167             wall_u(j,i) = 1.0_wp   ! south wall (location of adjacent fluid)
1168             fym(j,i)    = 1.0_wp
1169             fyp(j,i)    = 0.0_wp
[1]1170          ENDIF
1171!
1172!--       v-component
1173          IF ( nzb_v_outer(j,i) > nzb_v_outer(j,i+1) )  THEN
[1353]1174             wall_v(j,i) = 1.0_wp   ! rigth wall (location of adjacent fluid)
1175             fxm(j,i)    = 0.0_wp
1176             fxp(j,i)    = 1.0_wp
[1]1177          ELSEIF ( nzb_v_outer(j,i) > nzb_v_outer(j,i-1) )  THEN
[1353]1178             wall_v(j,i) = 1.0_wp   ! left wall (location of adjacent fluid)
1179             fxm(j,i)    = 1.0_wp
1180             fxp(j,i)    = 0.0_wp
[1]1181          ENDIF
1182!
1183!--       w-component, also used for scalars, separate arrays for shear
1184!--       production of tke
1185          IF ( nzb_w_outer(j,i) > nzb_w_outer(j+1,i) )  THEN
[1353]1186             wall_e_y(j,i) =  1.0_wp   ! north wall (location of adjacent fluid)
1187             wall_w_y(j,i) =  1.0_wp
1188             fwym(j,i)     =  0.0_wp
1189             fwyp(j,i)     =  1.0_wp
[1]1190          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j-1,i) )  THEN
[1353]1191             wall_e_y(j,i) = -1.0_wp   ! south wall (location of adjacent fluid)
1192             wall_w_y(j,i) =  1.0_wp
1193             fwym(j,i)     =  1.0_wp
1194             fwyp(j,i)     =  0.0_wp
[1]1195          ENDIF
1196          IF ( nzb_w_outer(j,i) > nzb_w_outer(j,i+1) )  THEN
[1353]1197             wall_e_x(j,i) =  1.0_wp   ! right wall (location of adjacent fluid)
1198             wall_w_x(j,i) =  1.0_wp
1199             fwxm(j,i)     =  0.0_wp
1200             fwxp(j,i)     =  1.0_wp
[1]1201          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j,i-1) )  THEN
[1353]1202             wall_e_x(j,i) = -1.0_wp   ! left wall (location of adjacent fluid)
1203             wall_w_x(j,i) =  1.0_wp
1204             fwxm(j,i)     =  1.0_wp
1205             fwxp(j,i)     =  0.0_wp
[1]1206          ENDIF
1207!
1208!--       Wall and corner locations inside buildings for limitation of
1209!--       near-wall mixing length l_wall
1210          IF ( nzb_s_inner(j,i) > nzb_s_inner(j+1,i) )  THEN
1211
1212             wall_n(j,i) = nzb_s_inner(j+1,i) + 1            ! North wall
1213
1214             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1215                corner_nl(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northleft corner
1216                                      nzb_s_inner(j,i-1) ) + 1
1217             ENDIF
1218
1219             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1220                corner_nr(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northright corner
1221                                      nzb_s_inner(j,i+1) ) + 1
1222             ENDIF
1223
1224          ENDIF
1225
1226          IF ( nzb_s_inner(j,i) > nzb_s_inner(j-1,i) )  THEN
1227
1228             wall_s(j,i) = nzb_s_inner(j-1,i) + 1            ! South wall
1229             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1230                corner_sl(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southleft corner
1231                                      nzb_s_inner(j,i-1) ) + 1
1232             ENDIF
1233
1234             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1235                corner_sr(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southright corner
1236                                      nzb_s_inner(j,i+1) ) + 1
1237             ENDIF
1238
1239          ENDIF
1240
1241          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1242             wall_l(j,i) = nzb_s_inner(j,i-1) + 1            ! Left wall
1243          ENDIF
1244
1245          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1246             wall_r(j,i) = nzb_s_inner(j,i+1) + 1            ! Right wall
1247          ENDIF
1248
1249       ENDDO
1250    ENDDO
1251!
[1931]1252!-- Calculate wall flag arrays for the multigrid method.
1253!-- Please note, wall flags are only applied in the not cache-optimized
1254!-- version.
1255    IF ( psolver == 'multigrid_noopt' )  THEN
[1968]1256
[114]1257!
[1968]1258!--    Gridpoint increment of the current level.
[114]1259       inc = 1
1260       DO  l = maximum_grid_level, 1 , -1
[1968]1261!
1262!--       Set grid_level as it is required for exchange_horiz_2d_int
1263          grid_level = l
[114]1264
1265          nxl_l = nxl_mg(l)
1266          nxr_l = nxr_mg(l)
1267          nys_l = nys_mg(l)
1268          nyn_l = nyn_mg(l)
1269          nzt_l = nzt_mg(l)
1270!
1271!--       Assign the flag level to be calculated
1272          SELECT CASE ( l )
1273             CASE ( 1 )
1274                flags => wall_flags_1
1275             CASE ( 2 )
1276                flags => wall_flags_2
1277             CASE ( 3 )
1278                flags => wall_flags_3
1279             CASE ( 4 )
1280                flags => wall_flags_4
1281             CASE ( 5 )
1282                flags => wall_flags_5
1283             CASE ( 6 )
1284                flags => wall_flags_6
1285             CASE ( 7 )
1286                flags => wall_flags_7
1287             CASE ( 8 )
1288                flags => wall_flags_8
1289             CASE ( 9 )
1290                flags => wall_flags_9
1291             CASE ( 10 )
1292                flags => wall_flags_10
1293          END SELECT
1294
1295!
1296!--       Depending on the grid level, set the respective bits in case of
1297!--       neighbouring walls
1298!--       Bit 0:  wall to the bottom
1299!--       Bit 1:  wall to the top (not realized in remaining PALM code so far)
1300!--       Bit 2:  wall to the south
1301!--       Bit 3:  wall to the north
1302!--       Bit 4:  wall to the left
1303!--       Bit 5:  wall to the right
[116]1304!--       Bit 6:  inside building
[114]1305
1306          flags = 0
1307
[927]1308!
1309!--       In case of masking method, flags are not set and multigrid method
1310!--       works like FFT-solver
[1931]1311          IF ( .NOT. masking_method )  THEN
[927]1312
[1968]1313!
1314!--          Allocate temporary array for topography heights on coarser grid
1315!--          level. Please note, 2 ghoist points are required, in order to
1316!--          calculate flags() on the interior ghost point.
1317             ALLOCATE( nzb_tmp(nys_l-2:nyn_l+2,nxl_l-2:nxr_l+2) )
1318             nzb_tmp = 0
1319             
1320             DO  i = nxl_l, nxr_l
1321                DO  j = nys_l, nyn_l
1322                   nzb_tmp(j,i) = nzb_local(j*inc,i*inc)
1323                ENDDO
1324             ENDDO
1325!
1326!--          Exchange ghost points on respective multigrid level. 2 ghost points
1327!--          are required, in order to calculate flags on
1328!--          nys_l-1 / nyn_l+1 / nxl_l-1 / nxr_l+1. The alternative would be to
1329!--          exchange 3D-INTEGER array flags on the respective multigrid level.
1330             CALL exchange_horiz_2d_int( nzb_tmp, nys_l, nyn_l, nxl_l, nxr_l, 2 )
1331!
1332!--          Set non-cyclic boundary conditions on respective multigrid level
1333             IF ( .NOT. bc_ns_cyc )  THEN
1334                IF ( nys == 0  )  THEN
1335                   nzb_tmp(-2,:) = nzb_tmp(0,:)
1336                   nzb_tmp(-1,:) = nzb_tmp(0,:)
1337                ENDIF
1338                IF ( nyn == ny )  THEN
1339                   nzb_tmp(ny+2,:) = nzb_tmp(ny,:)
1340                   nzb_tmp(ny+1,:) = nzb_tmp(ny,:)
1341                ENDIF
1342             ENDIF
1343             IF ( .NOT. bc_lr_cyc )  THEN
1344                IF ( nxl == 0  )  THEN
1345                   nzb_tmp(:,-2) = nzb_tmp(:,0)
1346                   nzb_tmp(:,-1) = nzb_tmp(:,0)
1347                ENDIF
1348                IF ( nxr == nx )  THEN
1349                   nzb_tmp(:,nx+1) = nzb_tmp(:,nx)   
1350                   nzb_tmp(:,nx+2) = nzb_tmp(:,nx)     
1351                ENDIF       
1352             ENDIF
1353                       
[927]1354             DO  i = nxl_l-1, nxr_l+1
1355                DO  j = nys_l-1, nyn_l+1
[1968]1356                   DO  k = nzb, nzt_l+1     
[114]1357!
[927]1358!--                   Inside/outside building (inside building does not need
1359!--                   further tests for walls)
[1968]1360                      IF ( k*inc <= nzb_tmp(j,i) )  THEN
[114]1361
[927]1362                         flags(k,j,i) = IBSET( flags(k,j,i), 6 )
[114]1363
[927]1364                      ELSE
[114]1365!
[927]1366!--                      Bottom wall
[1968]1367                         IF ( (k-1)*inc <= nzb_tmp(j,i) )  THEN
[927]1368                            flags(k,j,i) = IBSET( flags(k,j,i), 0 )
1369                         ENDIF
[114]1370!
[927]1371!--                      South wall
[1968]1372                         IF ( k*inc <= nzb_tmp(j-1,i) )  THEN
[927]1373                            flags(k,j,i) = IBSET( flags(k,j,i), 2 )
1374                         ENDIF
[114]1375!
[927]1376!--                      North wall
[1968]1377                         IF ( k*inc <= nzb_tmp(j+1,i) )  THEN
[927]1378                            flags(k,j,i) = IBSET( flags(k,j,i), 3 )
1379                         ENDIF
[114]1380!
[927]1381!--                      Left wall
[1968]1382                         IF ( k*inc <= nzb_tmp(j,i-1) )  THEN
[927]1383                            flags(k,j,i) = IBSET( flags(k,j,i), 4 )
1384                         ENDIF
[114]1385!
[927]1386!--                      Right wall
[1968]1387                         IF ( k*inc <= nzb_tmp(j,i+1) )  THEN
[927]1388                            flags(k,j,i) = IBSET( flags(k,j,i), 5 )
1389                         ENDIF
1390
[114]1391                      ENDIF
1392                           
[927]1393                   ENDDO
[114]1394                ENDDO
1395             ENDDO
1396
[1968]1397             DEALLOCATE( nzb_tmp )
1398
[927]1399          ENDIF
1400
[114]1401          inc = inc * 2
1402
1403       ENDDO
[1968]1404!
1405!--    Reset grid_level to "normal" grid
1406       grid_level = 0
1407       
[114]1408    ENDIF
[861]1409!
[1942]1410!-- Allocate flags needed for masking walls. Even though these flags are only
1411!-- required in the ws-scheme, the arrays need to be allocated as they are
1412!-- used in OpenACC directives.
[1677]1413    ALLOCATE( wall_flags_0(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
1414              wall_flags_00(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1221]1415    wall_flags_0  = 0
1416    wall_flags_00 = 0
[114]1417!
[1942]1418!-- Init flags for ws-scheme to degrade order near walls
1419    IF ( momentum_advec == 'ws-scheme'  .OR.  scalar_advec == 'ws-scheme'  .OR.&
1420         scalar_advec   == 'ws-scheme-mono' )  THEN
1421       CALL ws_init_flags
[861]1422    ENDIF
1423
1424!
[1]1425!-- In case of topography: limit near-wall mixing length l_wall further:
1426!-- Go through all points of the subdomain one by one and look for the closest
1427!-- surface
1428    IF ( TRIM(topography) /= 'flat' )  THEN
1429       DO  i = nxl, nxr
1430          DO  j = nys, nyn
1431
1432             nzb_si = nzb_s_inner(j,i)
1433             vi     = vertical_influence(nzb_si)
1434
1435             IF ( wall_n(j,i) > 0 )  THEN
1436!
1437!--             North wall (y distance)
1438                DO  k = wall_n(j,i), nzb_si
[1353]1439                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i), 0.5_wp * dy )
[1]1440                ENDDO
1441!
1442!--             Above North wall (yz distance)
1443                DO  k = nzb_si + 1, nzb_si + vi
[1353]1444                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i),                     &
1445                                          SQRT( 0.25_wp * dy**2 +              &
[1]1446                                          ( zu(k) - zw(nzb_si) )**2 ) )
1447                ENDDO
1448!
1449!--             Northleft corner (xy distance)
1450                IF ( corner_nl(j,i) > 0 )  THEN
1451                   DO  k = corner_nl(j,i), nzb_si
1452                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1), &
[1353]1453                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
[1]1454                   ENDDO
1455!
1456!--                Above Northleft corner (xyz distance)
1457                   DO  k = nzb_si + 1, nzb_si + vi
[1353]1458                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1),              &
1459                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1460                                            ( zu(k) - zw(nzb_si) )**2 ) )
[1]1461                   ENDDO
1462                ENDIF
1463!
1464!--             Northright corner (xy distance)
1465                IF ( corner_nr(j,i) > 0 )  THEN
1466                   DO  k = corner_nr(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 northright 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             ENDIF
1479
1480             IF ( wall_s(j,i) > 0 )  THEN
1481!
1482!--             South wall (y distance)
1483                DO  k = wall_s(j,i), nzb_si
[1353]1484                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i), 0.5_wp * dy )
[1]1485                ENDDO
1486!
1487!--             Above south wall (yz distance)
[1353]1488                DO  k = nzb_si + 1, nzb_si + vi
1489                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i),                     &
1490                                          SQRT( 0.25_wp * dy**2 +              &
[1]1491                                          ( zu(k) - zw(nzb_si) )**2 ) )
1492                ENDDO
1493!
1494!--             Southleft corner (xy distance)
1495                IF ( corner_sl(j,i) > 0 )  THEN
1496                   DO  k = corner_sl(j,i), nzb_si
[1353]1497                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1498                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
[1]1499                   ENDDO
1500!
1501!--                Above southleft corner (xyz distance)
1502                   DO  k = nzb_si + 1, nzb_si + vi
[1353]1503                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1504                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1505                                            ( zu(k) - zw(nzb_si) )**2 ) )
[1]1506                   ENDDO
1507                ENDIF
1508!
1509!--             Southright corner (xy distance)
1510                IF ( corner_sr(j,i) > 0 )  THEN
1511                   DO  k = corner_sr(j,i), nzb_si
[1353]1512                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1513                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
[1]1514                   ENDDO
1515!
1516!--                Above southright corner (xyz distance)
1517                   DO  k = nzb_si + 1, nzb_si + vi
[1353]1518                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1519                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1520                                            ( zu(k) - zw(nzb_si) )**2 ) )
[1]1521                   ENDDO
1522                ENDIF
1523
1524             ENDIF
1525
1526             IF ( wall_l(j,i) > 0 )  THEN
1527!
1528!--             Left wall (x distance)
1529                DO  k = wall_l(j,i), nzb_si
[1353]1530                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1), 0.5_wp * dx )
[1]1531                ENDDO
1532!
1533!--             Above left wall (xz distance)
1534                DO  k = nzb_si + 1, nzb_si + vi
[1353]1535                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1),                     &
1536                                       SQRT( 0.25_wp * dx**2 +                 &
1537                                       ( zu(k) - zw(nzb_si) )**2 ) )
[1]1538                ENDDO
1539             ENDIF
1540
1541             IF ( wall_r(j,i) > 0 )  THEN
1542!
1543!--             Right wall (x distance)
1544                DO  k = wall_r(j,i), nzb_si
[1353]1545                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1), 0.5_wp * dx )
[1]1546                ENDDO
1547!
1548!--             Above right wall (xz distance)
1549                DO  k = nzb_si + 1, nzb_si + vi
[1353]1550                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1),                     &
1551                                          SQRT( 0.25_wp * dx**2 +              &
[1]1552                                          ( zu(k) - zw(nzb_si) )**2 ) )
1553                ENDDO
1554
1555             ENDIF
1556
1557          ENDDO
1558       ENDDO
1559
1560    ENDIF
1561
1562!
1563!-- Multiplication with wall_adjustment_factor
1564    l_wall = wall_adjustment_factor * l_wall
1565
1566!
[709]1567!-- Set lateral boundary conditions for l_wall
[667]1568    CALL exchange_horiz( l_wall, nbgp )
1569
[1]1570    DEALLOCATE( corner_nl, corner_nr, corner_sl, corner_sr, nzb_local, &
[1968]1571                vertical_influence, wall_l, wall_n, wall_r, wall_s )
[1]1572
1573
1574 END SUBROUTINE init_grid
Note: See TracBrowser for help on using the repository browser.