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

Last change on this file since 3049 was 3049, checked in by Giersch, 6 years ago

Revision history corrected

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