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

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

changes from last commit documented

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