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

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

last commit documented

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