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

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

Rename multigrid into multigrid_noopt and multigrid_fast into multigrid, subroutines poismg is renamed into poismg_noopt and poismg_fast_mod into poismg_mod

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