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

Last change on this file since 1743 was 1743, checked in by raasch, 8 years ago

Bugfix for calculation of nzb_s_outer and nzb_u_outer at north boundary of total domain

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