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

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

last commit documented

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