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

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

get topograpyhy top index via function call

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