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