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

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

last commit documented

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