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

Last change on this file since 22 was 19, checked in by raasch, 18 years ago

preliminary version of modified boundary conditions at top

  • Property svn:keywords set to Id
File size: 29.0 KB
Line 
1 SUBROUTINE init_grid
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6! Setting of nzt_diff
7!
8! Former revisions:
9! -----------------
10! $Id: init_grid.f90 19 2007-02-23 04:53:48Z raasch $
11! RCS Log replace by Id keyword, revision history cleaned up
12!
13! Revision 1.17  2006/08/22 14:00:05  raasch
14! +dz_max to limit vertical stretching,
15! bugfix in index array initialization for line- or point-like topography
16! structures
17!
18! Revision 1.1  1997/08/11 06:17:45  raasch
19! Initial revision (Testversion)
20!
21!
22! Description:
23! ------------
24! Creating grid depending constants
25!------------------------------------------------------------------------------!
26
27    USE arrays_3d
28    USE control_parameters
29    USE grid_variables
30    USE indices
31    USE pegrid
32
33    IMPLICIT NONE
34
35    INTEGER ::  bh, blx, bly, bxl, bxr, byn, bys, i, i_center, j, j_center, k, &
36                l, nzb_si, nzt_l, vi
37
38    INTEGER, DIMENSION(:), ALLOCATABLE   ::  vertical_influence
39
40    INTEGER, DIMENSION(:,:), ALLOCATABLE ::  corner_nl, corner_nr, corner_sl,  &
41                                             corner_sr, wall_l, wall_n, wall_r,&
42                                             wall_s, nzb_local, nzb_tmp
43
44    REAL    ::  dx_l, dy_l, dz_stretched
45
46    REAL, DIMENSION(0:ny,0:nx)          ::  topo_height
47
48    REAL, DIMENSION(:,:,:), ALLOCATABLE ::  distance
49
50!
51!-- Allocate grid arrays
52    ALLOCATE( ddzu(1:nzt+1), ddzw(1:nzt+1), dd2zu(1:nzt), dzu(1:nzt+1), &
53              dzw(1:nzt+1), l_grid(1:nzt), zu(0:nzt+1), zw(0:nzt+1) )
54
55!
56!-- Compute height of u-levels from constant grid length and dz stretch factors
57    IF ( dz == -1.0 )  THEN
58       IF ( myid == 0 )  PRINT*,'+++ init_grid:  missing dz'
59       CALL local_stop
60    ELSEIF ( dz <= 0.0 )  THEN
61       IF ( myid == 0 )  PRINT*,'+++ init_grid:  dz=',dz,' <= 0.0'
62       CALL local_stop
63    ENDIF
64!
65!-- Since the w-level lies on the surface, the first u-level (staggered!) lies
66!-- below the surface (used for "mirror" boundary condition).
67!-- The first u-level above the surface corresponds to the top of the
68!-- Prandtl-layer.
69    zu(0) = - dz * 0.5
70    zu(1) =   dz * 0.5
71
72    dz_stretch_level_index = nzt+1
73    dz_stretched = dz
74    DO  k = 2, nzt+1
75       IF ( dz_stretch_level <= zu(k-1)  .AND.  dz_stretched < dz_max )  THEN
76          dz_stretched = dz_stretched * dz_stretch_factor
77          dz_stretched = MIN( dz_stretched, dz_max )
78          IF ( dz_stretch_level_index == nzt+1 )  dz_stretch_level_index = k-1
79       ENDIF
80       zu(k) = zu(k-1) + dz_stretched
81    ENDDO
82
83!
84!-- Compute the u-levels. They are always staggered half-way between the
85!-- corresponding w-levels. The top w-level is extrapolated linearly.
86    zw(0) = 0.0
87    DO  k = 1, nzt
88       zw(k) = ( zu(k) + zu(k+1) ) * 0.5
89    ENDDO
90    zw(nzt+1) = zw(nzt) + 2.0 * ( zu(nzt+1) - zw(nzt) )
91
92!
93!-- Compute grid lengths.
94    DO  k = 1, nzt+1
95       dzu(k)  = zu(k) - zu(k-1)
96       ddzu(k) = 1.0 / dzu(k)
97       dzw(k)  = zw(k) - zw(k-1)
98       ddzw(k) = 1.0 / dzw(k)
99    ENDDO
100
101    DO  k = 1, nzt
102       dd2zu(k) = 1.0 / ( dzu(k) + dzu(k+1) )
103    ENDDO
104
105!
106!-- In case of multigrid method, compute grid lengths and grid factors for the
107!-- grid levels
108    IF ( psolver == 'multigrid' )  THEN
109
110       ALLOCATE( ddx2_mg(maximum_grid_level), ddy2_mg(maximum_grid_level), &
111                 dzu_mg(nzb+1:nzt+1,maximum_grid_level),                   &
112                 dzw_mg(nzb+1:nzt+1,maximum_grid_level),                   &
113                 f1_mg(nzb+1:nzt,maximum_grid_level),                      &
114                 f2_mg(nzb+1:nzt,maximum_grid_level),                      &
115                 f3_mg(nzb+1:nzt,maximum_grid_level) )
116
117       dzu_mg(:,maximum_grid_level) = dzu
118       dzw_mg(:,maximum_grid_level) = dzw
119       nzt_l = nzt
120       DO  l = maximum_grid_level-1, 1, -1
121           dzu_mg(nzb+1,l) = 2.0 * dzu_mg(nzb+1,l+1)
122           dzw_mg(nzb+1,l) = 2.0 * dzw_mg(nzb+1,l+1)
123           nzt_l = nzt_l / 2
124           DO  k = 2, nzt_l+1
125              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
126              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
127           ENDDO
128       ENDDO
129
130       nzt_l = nzt
131       dx_l  = dx
132       dy_l  = dy
133       DO  l = maximum_grid_level, 1, -1
134          ddx2_mg(l) = 1.0 / dx_l**2
135          ddy2_mg(l) = 1.0 / dy_l**2
136          DO  k = nzb+1, nzt_l
137             f2_mg(k,l) = 1.0 / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
138             f3_mg(k,l) = 1.0 / ( dzu_mg(k,l)   * dzw_mg(k,l) )
139             f1_mg(k,l) = 2.0 * ( ddx2_mg(l) + ddy2_mg(l) ) + &
140                          f2_mg(k,l) + f3_mg(k,l)
141          ENDDO
142          nzt_l = nzt_l / 2
143          dx_l  = dx_l * 2.0
144          dy_l  = dy_l * 2.0
145       ENDDO
146
147    ENDIF
148
149!
150!-- Compute the reciprocal values of the horizontal grid lengths.
151    ddx = 1.0 / dx
152    ddy = 1.0 / dy
153    dx2 = dx * dx
154    dy2 = dy * dy
155    ddx2 = 1.0 / dx2
156    ddy2 = 1.0 / dy2
157
158!
159!-- Compute the grid-dependent mixing length.
160    DO  k = 1, nzt
161       l_grid(k)  = ( dx * dy * dzw(k) )**0.33333333333333
162    ENDDO
163
164!
165!-- Allocate outer and inner index arrays for topography and set
166!-- defaults
167    ALLOCATE( corner_nl(nys:nyn,nxl:nxr), corner_nr(nys:nyn,nxl:nxr), &
168              corner_sl(nys:nyn,nxl:nxr), corner_sr(nys:nyn,nxl:nxr), &
169              nzb_local(-1:ny+1,-1:nx+1), nzb_tmp(-1:ny+1,-1:nx+1),   &
170              wall_l(nys:nyn,nxl:nxr), wall_n(nys:nyn,nxl:nxr),       &
171              wall_r(nys:nyn,nxl:nxr), wall_s(nys:nyn,nxl:nxr) )
172    ALLOCATE( fwxm(nys-1:nyn+1,nxl-1:nxr+1), fwxp(nys-1:nyn+1,nxl-1:nxr+1), &
173              fwym(nys-1:nyn+1,nxl-1:nxr+1), fwyp(nys-1:nyn+1,nxl-1:nxr+1), &
174              fxm(nys-1:nyn+1,nxl-1:nxr+1), fxp(nys-1:nyn+1,nxl-1:nxr+1),   &
175              fym(nys-1:nyn+1,nxl-1:nxr+1), fyp(nys-1:nyn+1,nxl-1:nxr+1),   &
176              nzb_s_inner(nys-1:nyn+1,nxl-1:nxr+1),                         &
177              nzb_s_outer(nys-1:nyn+1,nxl-1:nxr+1),                         &
178              nzb_u_inner(nys-1:nyn+1,nxl-1:nxr+1),                         &
179              nzb_u_outer(nys-1:nyn+1,nxl-1:nxr+1),                         &
180              nzb_v_inner(nys-1:nyn+1,nxl-1:nxr+1),                         &
181              nzb_v_outer(nys-1:nyn+1,nxl-1:nxr+1),                         &
182              nzb_w_inner(nys-1:nyn+1,nxl-1:nxr+1),                         &
183              nzb_w_outer(nys-1:nyn+1,nxl-1:nxr+1),                         &
184              nzb_diff_s_inner(nys-1:nyn+1,nxl-1:nxr+1),                    &
185              nzb_diff_s_outer(nys-1:nyn+1,nxl-1:nxr+1),                    &
186              nzb_diff_u(nys-1:nyn+1,nxl-1:nxr+1),                          &
187              nzb_diff_v(nys-1:nyn+1,nxl-1:nxr+1),                          &
188              nzb_2d(nys-1:nyn+1,nxl-1:nxr+1),                              &
189              wall_e_x(nys-1:nyn+1,nxl-1:nxr+1),                            &
190              wall_e_y(nys-1:nyn+1,nxl-1:nxr+1),                            &
191              wall_u(nys-1:nyn+1,nxl-1:nxr+1),                              &
192              wall_v(nys-1:nyn+1,nxl-1:nxr+1),                              &
193              wall_w_x(nys-1:nyn+1,nxl-1:nxr+1),                            &
194              wall_w_y(nys-1:nyn+1,nxl-1:nxr+1) )
195
196    ALLOCATE( l_wall(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
197
198    nzb_s_inner = nzb;  nzb_s_outer = nzb
199    nzb_u_inner = nzb;  nzb_u_outer = nzb
200    nzb_v_inner = nzb;  nzb_v_outer = nzb
201    nzb_w_inner = nzb;  nzb_w_outer = nzb
202
203!
204!-- Define vertical gridpoint from (or to) which on the usual finite difference
205!-- form (which does not use surface fluxes) is applied
206    IF ( prandtl_layer  .OR.  use_surface_fluxes )  THEN
207       nzb_diff = nzb + 2
208    ELSE
209       nzb_diff = nzb + 1
210    ENDIF
211    IF ( use_top_fluxes )  THEN
212       nzt_diff = nzt - 1
213    ELSE
214       nzt_diff = nzt
215    ENDIF
216
217    nzb_diff_s_inner = nzb_diff;  nzb_diff_s_outer = nzb_diff
218    nzb_diff_u = nzb_diff;  nzb_diff_v = nzb_diff
219
220    wall_e_x = 0.0;  wall_e_y = 0.0;  wall_u = 0.0;  wall_v = 0.0
221    wall_w_x = 0.0;  wall_w_y = 0.0
222    fwxp = 1.0;  fwxm = 1.0;  fwyp = 1.0;  fwym = 1.0
223    fxp  = 1.0;  fxm  = 1.0;  fyp  = 1.0;  fym  = 1.0
224
225!
226!-- Initialize near-wall mixing length l_wall only in the vertical direction
227!-- for the moment,
228!-- multiplication with wall_adjustment_factor near the end of this routine
229    l_wall(nzb,:,:)   = l_grid(1)
230    DO  k = nzb+1, nzt
231       l_wall(k,:,:)  = l_grid(k)
232    ENDDO
233    l_wall(nzt+1,:,:) = l_grid(nzt)
234
235    ALLOCATE ( vertical_influence(nzb:nzt) )
236    DO  k = 1, nzt
237       vertical_influence(k) = MIN ( INT( l_grid(k) / &
238                     ( wall_adjustment_factor * dzw(k) ) + 0.5 ), nzt - k )
239    ENDDO
240
241    DO  k = 1, MAXVAL( nzb_s_inner )
242       IF ( l_grid(k) > 1.5 * dx * wall_adjustment_factor .OR.  &
243            l_grid(k) > 1.5 * dy * wall_adjustment_factor )  THEN
244          IF ( myid == 0 )  THEN
245             PRINT*, '+++ WARNING: grid anisotropy exceeds '// &
246                           'threshold given by only local'
247             PRINT*, '             horizontal reduction of near_wall '// &
248                           'mixing length l_wall'
249             PRINT*, '             starting from height level k = ', k, '.'
250          ENDIF
251          EXIT
252       ENDIF
253    ENDDO
254    vertical_influence(0) = vertical_influence(1)
255
256    DO  i = nxl-1, nxr+1
257       DO  j = nys-1, nyn+1
258          DO  k = nzb_s_inner(j,i) + 1, &
259                  nzb_s_inner(j,i) + vertical_influence(nzb_s_inner(j,i))
260             l_wall(k,j,i) = zu(k) - zw(nzb_s_inner(j,i))
261          ENDDO
262       ENDDO
263    ENDDO
264
265!
266!-- Set outer and inner index arrays for non-flat topography.
267!-- Here consistency checks concerning domain size and periodicity are
268!-- necessary.
269!-- Within this SELECT CASE structure only nzb_local is initialized
270!-- individually depending on the chosen topography type, all other index
271!-- arrays are initialized further below.
272    SELECT CASE ( TRIM( topography ) )
273
274       CASE ( 'flat' )
275!
276!--       No actions necessary
277
278       CASE ( 'single_building' )
279!
280!--       Single rectangular building, by default centered in the middle of the
281!--       total domain
282          blx = NINT( building_length_x / dx )
283          bly = NINT( building_length_y / dy )
284          bh  = NINT( building_height / dz )
285
286          IF ( building_wall_left == 9999999.9 )  THEN
287             building_wall_left = ( nx + 1 - blx ) / 2 * dx
288          ENDIF
289          bxl = NINT( building_wall_left / dx )
290          bxr = bxl + blx
291
292          IF ( building_wall_south == 9999999.9 )  THEN
293             building_wall_south = ( ny + 1 - bly ) / 2 * dy
294          ENDIF
295          bys = NINT( building_wall_south / dy )
296          byn = bys + bly
297
298!
299!--       Building size has to meet some requirements
300          IF ( ( bxl < 1 ) .OR. ( bxr > nx-1 ) .OR. ( bxr < bxl+3 ) .OR.  &
301               ( bys < 1 ) .OR. ( byn > ny-1 ) .OR. ( byn < bys+3 ) )  THEN
302             IF ( myid == 0 )  THEN
303                PRINT*, '+++ init_grid: inconsistent building parameters:'
304                PRINT*, '               bxl=', bxl, 'bxr=', bxr, 'bys=', bys, &
305                                      'byn=', byn, 'nx=', nx, 'ny=', ny
306             ENDIF
307             CALL local_stop
308          ENDIF
309
310!
311!--       Set the individual index arrays for all velocity components and
312!--       scalars, taking into account the staggered grid. The horizontal
313!--       wind component normal to a wall defines the position of the wall, and
314!--       in the respective direction the building is as long as specified in
315!--       building_length_?, but in the other horizontal direction (for w and s
316!--       in both horizontal directions) the building appears shortened by one
317!--       grid length due to the staggered grid.
318          nzb_local = 0
319          nzb_local(bys:byn-1,bxl:bxr-1) = bh
320
321       CASE ( 'read_from_file' )
322!
323!--       Arbitrary irregular topography data in PALM format (exactly matching
324!--       the grid size and total domain size)
325          OPEN( 90, FILE='TOPOGRAPHY_DATA', STATUS='OLD', FORM='FORMATTED',  &
326               ERR=10 )
327          DO  j = ny, 0, -1
328             READ( 90, *, ERR=11, END=11 )  ( topo_height(j,i), i = 0, nx )
329          ENDDO
330!
331!--       Calculate the index height of the topography
332          DO  i = 0, nx
333             DO  j = 0, ny
334                nzb_local(j,i) = NINT( topo_height(j,i) / dz )
335             ENDDO
336          ENDDO
337          nzb_local(-1,0:nx)   = nzb_local(ny,0:nx)
338          nzb_local(ny+1,0:nx) = nzb_local(0,0:nx)
339          nzb_local(:,-1)      = nzb_local(:,nx)
340          nzb_local(:,nx+1)    = nzb_local(:,0)
341     
342          GOTO 12
343
344 10       IF ( myid == 0 )  THEN
345             PRINT*, '+++ init_grid: file TOPOGRAPHY_DATA does not exist'
346          ENDIF
347          CALL local_stop
348
349 11       IF ( myid == 0 )  THEN
350             PRINT*, '+++ init_grid: errors in file TOPOGRAPHY_DATA'
351          ENDIF
352          CALL local_stop
353
354 12       CLOSE( 90 )
355
356       CASE DEFAULT
357!
358!--       The DEFAULT case is reached either if the parameter topography
359!--       contains a wrong character string or if the user has coded a special
360!--       case in the user interface. There, the subroutine user_init_grid
361!--       checks which of these two conditions applies.
362          CALL user_init_grid( nzb_local )
363
364    END SELECT
365
366!
367!-- Consistency checks and index array initialization are only required for
368!-- non-flat topography
369    IF ( TRIM( topography ) /= 'flat' )  THEN
370
371!
372!--    Consistency checks
373       IF ( MINVAL( nzb_local ) < 0  .OR.  MAXVAL( nzb_local ) > nz + 1 )  THEN
374          IF ( myid == 0 )  THEN
375             PRINT*, '+++ init_grid: nzb_local values are outside the', &
376                          'model domain'
377             PRINT*, '               MINVAL( nzb_local ) = ', MINVAL(nzb_local)
378             PRINT*, '               MAXVAL( nzb_local ) = ', MAXVAL(nzb_local)
379          ENDIF
380          CALL local_stop
381       ENDIF
382
383       IF ( bc_lr == 'cyclic' )  THEN
384          IF ( ANY( nzb_local(:,-1) /= nzb_local(:,nx)   )  .OR. &
385               ANY( nzb_local(:,0)  /= nzb_local(:,nx+1) ) )  THEN
386             IF ( myid == 0 )  THEN
387                PRINT*, '+++ init_grid: nzb_local does not fulfill cyclic', &
388                        '               boundary condition in x-direction'
389             ENDIF
390             CALL local_stop
391          ENDIF
392       ENDIF
393       IF ( bc_ns == 'cyclic' )  THEN
394          IF ( ANY( nzb_local(-1,:) /= nzb_local(ny,:)   )  .OR. &
395               ANY( nzb_local(0,:)  /= nzb_local(ny+1,:) ) )  THEN
396             IF ( myid == 0 )  THEN
397                PRINT*, '+++ init_grid: nzb_local does not fulfill cyclic', &
398                        '               boundary condition in y-direction'
399             ENDIF
400             CALL local_stop
401          ENDIF
402       ENDIF
403
404!
405!--    Initialize index arrays nzb_s_inner and nzb_w_inner
406       nzb_s_inner = nzb_local(nys-1:nyn+1,nxl-1:nxr+1)
407       nzb_w_inner = nzb_local(nys-1:nyn+1,nxl-1:nxr+1)
408
409!
410!--    Initialize remaining index arrays:
411!--    first pre-initialize them with nzb_s_inner...
412       nzb_u_inner = nzb_s_inner
413       nzb_u_outer = nzb_s_inner
414       nzb_v_inner = nzb_s_inner
415       nzb_v_outer = nzb_s_inner
416       nzb_w_outer = nzb_s_inner
417       nzb_s_outer = nzb_s_inner
418
419!
420!--    ...then extend pre-initialized arrays in their according directions
421!--    based on nzb_local using nzb_tmp as a temporary global index array
422
423!
424!--    nzb_s_outer:
425!--    extend nzb_local east-/westwards first, then north-/southwards
426       nzb_tmp = nzb_local
427       DO  j = -1, ny + 1
428          DO  i = 0, nx
429             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i), &
430                                 nzb_local(j,i+1) )
431          ENDDO
432       ENDDO
433       DO  i = nxl, nxr
434          DO  j = nys, nyn
435             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i), &
436                                     nzb_tmp(j+1,i) )
437          ENDDO
438!
439!--       non-cyclic boundary conditions (overwritten by call of
440!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
441          IF ( nys == 0 )  THEN
442             j = -1
443             nzb_s_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
444          ENDIF
445          IF ( nys == ny )  THEN
446             j = ny + 1
447             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
448          ENDIF
449       ENDDO
450!
451!--    nzb_w_outer:
452!--    identical to nzb_s_outer
453       nzb_w_outer = nzb_s_outer
454
455!
456!--    nzb_u_inner:
457!--    extend nzb_local rightwards only
458       nzb_tmp = nzb_local
459       DO  j = -1, ny + 1
460          DO  i = 0, nx + 1
461             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i) )
462          ENDDO
463       ENDDO
464       nzb_u_inner = nzb_tmp(nys-1:nyn+1,nxl-1:nxr+1)
465
466!
467!--    nzb_u_outer:
468!--    extend current nzb_tmp (nzb_u_inner) north-/southwards
469       DO  i = nxl, nxr
470          DO  j = nys, nyn
471             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i), &
472                                     nzb_tmp(j+1,i) )
473          ENDDO
474!
475!--       non-cyclic boundary conditions (overwritten by call of
476!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
477          IF ( nys == 0 )  THEN
478             j = -1
479             nzb_u_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
480          ENDIF
481          IF ( nys == ny )  THEN
482             j = ny + 1
483             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
484          ENDIF
485       ENDDO
486
487!
488!--    nzb_v_inner:
489!--    extend nzb_local northwards only
490       nzb_tmp = nzb_local
491       DO  i = -1, nx + 1
492          DO  j = 0, ny + 1
493             nzb_tmp(j,i) = MAX( nzb_local(j-1,i), nzb_local(j,i) )
494          ENDDO
495       ENDDO
496       nzb_v_inner = nzb_tmp(nys-1:nyn+1,nxl-1:nxr+1)
497
498!
499!--    nzb_v_outer:
500!--    extend current nzb_tmp (nzb_v_inner) right-/leftwards
501       DO  j = nys, nyn
502          DO  i = nxl, nxr
503             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i), &
504                                     nzb_tmp(j,i+1) )
505          ENDDO
506!
507!--       non-cyclic boundary conditions (overwritten by call of
508!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
509          IF ( nxl == 0 )  THEN
510             i = -1
511             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i+1), nzb_tmp(j,i) )
512          ENDIF
513          IF ( nxr == nx )  THEN
514             i = nx + 1
515             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i) )
516          ENDIF
517       ENDDO
518
519!
520!--    Exchange of lateral boundary values (parallel computers) and cyclic
521!--    boundary conditions, if applicable.
522!--    Since nzb_s_inner and nzb_w_inner are derived directly from nzb_local
523!--    they do not require exchange and are not included here.
524       CALL exchange_horiz_2d_int( nzb_u_inner )
525       CALL exchange_horiz_2d_int( nzb_u_outer )
526       CALL exchange_horiz_2d_int( nzb_v_inner )
527       CALL exchange_horiz_2d_int( nzb_v_outer )
528       CALL exchange_horiz_2d_int( nzb_w_outer )
529       CALL exchange_horiz_2d_int( nzb_s_outer )
530
531    ENDIF
532
533!
534!-- Preliminary: to be removed after completion of the topography code!
535!-- Set the former default k index arrays nzb_2d
536    nzb_2d      = nzb
537
538!
539!-- Set the individual index arrays which define the k index from which on
540!-- the usual finite difference form (which does not use surface fluxes) is
541!-- applied
542    IF ( prandtl_layer  .OR.  use_surface_fluxes )  THEN
543       nzb_diff_u         = nzb_u_inner + 2
544       nzb_diff_v         = nzb_v_inner + 2
545       nzb_diff_s_inner   = nzb_s_inner + 2
546       nzb_diff_s_outer   = nzb_s_outer + 2
547    ELSE
548       nzb_diff_u         = nzb_u_inner + 1
549       nzb_diff_v         = nzb_v_inner + 1
550       nzb_diff_s_inner   = nzb_s_inner + 1
551       nzb_diff_s_outer   = nzb_s_outer + 1
552    ENDIF
553
554!
555!-- Calculation of wall switches and factors required by diffusion_u/v.f90 and
556!-- for limitation of near-wall mixing length l_wall further below
557    corner_nl = 0
558    corner_nr = 0
559    corner_sl = 0
560    corner_sr = 0
561    wall_l    = 0
562    wall_n    = 0
563    wall_r    = 0
564    wall_s    = 0
565
566    DO  i = nxl, nxr
567       DO  j = nys, nyn
568!
569!--       u-component
570          IF ( nzb_u_outer(j,i) > nzb_u_outer(j+1,i) )  THEN
571             wall_u(j,i) = 1.0   ! north wall (location of adjacent fluid)
572             fym(j,i)    = 0.0
573             fyp(j,i)    = 1.0
574          ELSEIF ( nzb_u_outer(j,i) > nzb_u_outer(j-1,i) )  THEN
575             wall_u(j,i) = 1.0   ! south wall (location of adjacent fluid)
576             fym(j,i)    = 1.0
577             fyp(j,i)    = 0.0
578          ENDIF
579!
580!--       v-component
581          IF ( nzb_v_outer(j,i) > nzb_v_outer(j,i+1) )  THEN
582             wall_v(j,i) = 1.0   ! rigth wall (location of adjacent fluid)
583             fxm(j,i)    = 0.0
584             fxp(j,i)    = 1.0
585          ELSEIF ( nzb_v_outer(j,i) > nzb_v_outer(j,i-1) )  THEN
586             wall_v(j,i) = 1.0   ! left wall (location of adjacent fluid)
587             fxm(j,i)    = 1.0
588             fxp(j,i)    = 0.0
589          ENDIF
590!
591!--       w-component, also used for scalars, separate arrays for shear
592!--       production of tke
593          IF ( nzb_w_outer(j,i) > nzb_w_outer(j+1,i) )  THEN
594             wall_e_y(j,i) =  1.0   ! north wall (location of adjacent fluid)
595             wall_w_y(j,i) =  1.0
596             fwym(j,i)     =  0.0
597             fwyp(j,i)     =  1.0
598          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j-1,i) )  THEN
599             wall_e_y(j,i) = -1.0   ! south wall (location of adjacent fluid)
600             wall_w_y(j,i) =  1.0
601             fwym(j,i)     =  1.0
602             fwyp(j,i)     =  0.0
603          ENDIF
604          IF ( nzb_w_outer(j,i) > nzb_w_outer(j,i+1) )  THEN
605             wall_e_x(j,i) =  1.0   ! right wall (location of adjacent fluid)
606             wall_w_x(j,i) =  1.0
607             fwxm(j,i)     =  0.0
608             fwxp(j,i)     =  1.0
609          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j,i-1) )  THEN
610             wall_e_x(j,i) = -1.0   ! left wall (location of adjacent fluid)
611             wall_w_x(j,i) =  1.0
612             fwxm(j,i)     =  1.0
613             fwxp(j,i)     =  0.0
614          ENDIF
615!
616!--       Wall and corner locations inside buildings for limitation of
617!--       near-wall mixing length l_wall
618          IF ( nzb_s_inner(j,i) > nzb_s_inner(j+1,i) )  THEN
619
620             wall_n(j,i) = nzb_s_inner(j+1,i) + 1            ! North wall
621
622             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
623                corner_nl(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northleft corner
624                                      nzb_s_inner(j,i-1) ) + 1
625             ENDIF
626
627             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
628                corner_nr(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northright corner
629                                      nzb_s_inner(j,i+1) ) + 1
630             ENDIF
631
632          ENDIF
633
634          IF ( nzb_s_inner(j,i) > nzb_s_inner(j-1,i) )  THEN
635
636             wall_s(j,i) = nzb_s_inner(j-1,i) + 1            ! South wall
637             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
638                corner_sl(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southleft corner
639                                      nzb_s_inner(j,i-1) ) + 1
640             ENDIF
641
642             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
643                corner_sr(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southright corner
644                                      nzb_s_inner(j,i+1) ) + 1
645             ENDIF
646
647          ENDIF
648
649          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
650             wall_l(j,i) = nzb_s_inner(j,i-1) + 1            ! Left wall
651          ENDIF
652
653          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
654             wall_r(j,i) = nzb_s_inner(j,i+1) + 1            ! Right wall
655          ENDIF
656
657       ENDDO
658    ENDDO
659
660!
661!-- In case of topography: limit near-wall mixing length l_wall further:
662!-- Go through all points of the subdomain one by one and look for the closest
663!-- surface
664    IF ( TRIM(topography) /= 'flat' )  THEN
665       DO  i = nxl, nxr
666          DO  j = nys, nyn
667
668             nzb_si = nzb_s_inner(j,i)
669             vi     = vertical_influence(nzb_si)
670
671             IF ( wall_n(j,i) > 0 )  THEN
672!
673!--             North wall (y distance)
674                DO  k = wall_n(j,i), nzb_si
675                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i), 0.5 * dy )
676                ENDDO
677!
678!--             Above North wall (yz distance)
679                DO  k = nzb_si + 1, nzb_si + vi
680                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i),     &
681                                          SQRT( 0.25 * dy**2 + &
682                                          ( zu(k) - zw(nzb_si) )**2 ) )
683                ENDDO
684!
685!--             Northleft corner (xy distance)
686                IF ( corner_nl(j,i) > 0 )  THEN
687                   DO  k = corner_nl(j,i), nzb_si
688                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1), &
689                                               0.5 * SQRT( dx**2 + dy**2 ) )
690                   ENDDO
691!
692!--                Above Northleft corner (xyz distance)
693                   DO  k = nzb_si + 1, nzb_si + vi
694                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1),             &
695                                               SQRT( 0.25 * (dx**2 + dy**2) + &
696                                               ( zu(k) - zw(nzb_si) )**2 ) )
697                   ENDDO
698                ENDIF
699!
700!--             Northright corner (xy distance)
701                IF ( corner_nr(j,i) > 0 )  THEN
702                   DO  k = corner_nr(j,i), nzb_si
703                       l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1), &
704                                                0.5 * SQRT( dx**2 + dy**2 ) )
705                   ENDDO
706!
707!--                Above northright corner (xyz distance)
708                   DO  k = nzb_si + 1, nzb_si + vi
709                      l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1), &
710                                               SQRT( 0.25 * (dx**2 + dy**2) + &
711                                               ( zu(k) - zw(nzb_si) )**2 ) )
712                   ENDDO
713                ENDIF
714             ENDIF
715
716             IF ( wall_s(j,i) > 0 )  THEN
717!
718!--             South wall (y distance)
719                DO  k = wall_s(j,i), nzb_si
720                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i), 0.5 * dy )
721                ENDDO
722!
723!--             Above south wall (yz distance)
724                DO  k = nzb_si + 1, &
725                        nzb_si + vi
726                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i),     &
727                                          SQRT( 0.25 * dy**2 + &
728                                          ( zu(k) - zw(nzb_si) )**2 ) )
729                ENDDO
730!
731!--             Southleft corner (xy distance)
732                IF ( corner_sl(j,i) > 0 )  THEN
733                   DO  k = corner_sl(j,i), nzb_si
734                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1), &
735                                               0.5 * SQRT( dx**2 + dy**2 ) )
736                   ENDDO
737!
738!--                Above southleft corner (xyz distance)
739                   DO  k = nzb_si + 1, nzb_si + vi
740                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),             &
741                                               SQRT( 0.25 * (dx**2 + dy**2) + &
742                                               ( zu(k) - zw(nzb_si) )**2 ) )
743                   ENDDO
744                ENDIF
745!
746!--             Southright corner (xy distance)
747                IF ( corner_sr(j,i) > 0 )  THEN
748                   DO  k = corner_sr(j,i), nzb_si
749                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1), &
750                                               0.5 * SQRT( dx**2 + dy**2 ) )
751                   ENDDO
752!
753!--                Above southright corner (xyz distance)
754                   DO  k = nzb_si + 1, nzb_si + vi
755                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),             &
756                                               SQRT( 0.25 * (dx**2 + dy**2) + &
757                                               ( zu(k) - zw(nzb_si) )**2 ) )
758                   ENDDO
759                ENDIF
760
761             ENDIF
762
763             IF ( wall_l(j,i) > 0 )  THEN
764!
765!--             Left wall (x distance)
766                DO  k = wall_l(j,i), nzb_si
767                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1), 0.5 * dx )
768                ENDDO
769!
770!--             Above left wall (xz distance)
771                DO  k = nzb_si + 1, nzb_si + vi
772                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1),     &
773                                          SQRT( 0.25 * dx**2 + &
774                                          ( zu(k) - zw(nzb_si) )**2 ) )
775                ENDDO
776             ENDIF
777
778             IF ( wall_r(j,i) > 0 )  THEN
779!
780!--             Right wall (x distance)
781                DO  k = wall_r(j,i), nzb_si
782                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1), 0.5 * dx )
783                ENDDO
784!
785!--             Above right wall (xz distance)
786                DO  k = nzb_si + 1, nzb_si + vi
787                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1),     &
788                                          SQRT( 0.25 * dx**2 + &
789                                          ( zu(k) - zw(nzb_si) )**2 ) )
790                ENDDO
791
792             ENDIF
793
794          ENDDO
795       ENDDO
796
797    ENDIF
798
799!
800!-- Multiplication with wall_adjustment_factor
801    l_wall = wall_adjustment_factor * l_wall
802
803!
804!-- Need to set lateral boundary conditions for l_wall
805    CALL exchange_horiz( l_wall, 0, 0 )
806
807    DEALLOCATE( corner_nl, corner_nr, corner_sl, corner_sr, nzb_local, &
808                nzb_tmp, vertical_influence, wall_l, wall_n, wall_r, wall_s )
809
810
811 END SUBROUTINE init_grid
Note: See TracBrowser for help on using the repository browser.