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

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

Bugfix, move setting of topography_grid_convention to init_grid, else, in case of restarts generic topography may not set properly if no value is prescribed

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