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

Last change on this file since 1763 was 1763, checked in by hellstea, 8 years ago

last commit documented

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