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

Last change on this file since 2022 was 2022, checked in by suehring, 7 years ago

last commit documented

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