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

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

Id keyword set as property for all *.f90 files

  • Property svn:keywords set to Id
File size: 28.8 KB
Line 
1 SUBROUTINE init_grid
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6!
7!
8! Former revisions:
9! -----------------
10! $Id: init_grid.f90 4 2007-02-13 11:33:16Z 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 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
212    nzb_diff_s_inner = nzb_diff;  nzb_diff_s_outer = nzb_diff
213    nzb_diff_u = nzb_diff;  nzb_diff_v = nzb_diff
214
215    wall_e_x = 0.0;  wall_e_y = 0.0;  wall_u = 0.0;  wall_v = 0.0
216    wall_w_x = 0.0;  wall_w_y = 0.0
217    fwxp = 1.0;  fwxm = 1.0;  fwyp = 1.0;  fwym = 1.0
218    fxp  = 1.0;  fxm  = 1.0;  fyp  = 1.0;  fym  = 1.0
219
220!
221!-- Initialize near-wall mixing length l_wall only in the vertical direction
222!-- for the moment,
223!-- multiplication with wall_adjustment_factor near the end of this routine
224    l_wall(nzb,:,:)   = l_grid(1)
225    DO  k = nzb+1, nzt
226       l_wall(k,:,:)  = l_grid(k)
227    ENDDO
228    l_wall(nzt+1,:,:) = l_grid(nzt)
229
230    ALLOCATE ( vertical_influence(nzb:nzt) )
231    DO  k = 1, nzt
232       vertical_influence(k) = MIN ( INT( l_grid(k) / &
233                     ( wall_adjustment_factor * dzw(k) ) + 0.5 ), nzt - k )
234    ENDDO
235
236    DO  k = 1, MAXVAL( nzb_s_inner )
237       IF ( l_grid(k) > 1.5 * dx * wall_adjustment_factor .OR.  &
238            l_grid(k) > 1.5 * dy * wall_adjustment_factor )  THEN
239          IF ( myid == 0 )  THEN
240             PRINT*, '+++ WARNING: grid anisotropy exceeds '// &
241                           'threshold given by only local'
242             PRINT*, '             horizontal reduction of near_wall '// &
243                           'mixing length l_wall'
244             PRINT*, '             starting from height level k = ', k, '.'
245          ENDIF
246          EXIT
247       ENDIF
248    ENDDO
249    vertical_influence(0) = vertical_influence(1)
250
251    DO  i = nxl-1, nxr+1
252       DO  j = nys-1, nyn+1
253          DO  k = nzb_s_inner(j,i) + 1, &
254                  nzb_s_inner(j,i) + vertical_influence(nzb_s_inner(j,i))
255             l_wall(k,j,i) = zu(k) - zw(nzb_s_inner(j,i))
256          ENDDO
257       ENDDO
258    ENDDO
259
260!
261!-- Set outer and inner index arrays for non-flat topography.
262!-- Here consistency checks concerning domain size and periodicity are
263!-- necessary.
264!-- Within this SELECT CASE structure only nzb_local is initialized
265!-- individually depending on the chosen topography type, all other index
266!-- arrays are initialized further below.
267    SELECT CASE ( TRIM( topography ) )
268
269       CASE ( 'flat' )
270!
271!--       No actions necessary
272
273       CASE ( 'single_building' )
274!
275!--       Single rectangular building, by default centered in the middle of the
276!--       total domain
277          blx = NINT( building_length_x / dx )
278          bly = NINT( building_length_y / dy )
279          bh  = NINT( building_height / dz )
280
281          IF ( building_wall_left == 9999999.9 )  THEN
282             building_wall_left = ( nx + 1 - blx ) / 2 * dx
283          ENDIF
284          bxl = NINT( building_wall_left / dx )
285          bxr = bxl + blx
286
287          IF ( building_wall_south == 9999999.9 )  THEN
288             building_wall_south = ( ny + 1 - bly ) / 2 * dy
289          ENDIF
290          bys = NINT( building_wall_south / dy )
291          byn = bys + bly
292
293!
294!--       Building size has to meet some requirements
295          IF ( ( bxl < 1 ) .OR. ( bxr > nx-1 ) .OR. ( bxr < bxl+3 ) .OR.  &
296               ( bys < 1 ) .OR. ( byn > ny-1 ) .OR. ( byn < bys+3 ) )  THEN
297             IF ( myid == 0 )  THEN
298                PRINT*, '+++ init_grid: inconsistent building parameters:'
299                PRINT*, '               bxl=', bxl, 'bxr=', bxr, 'bys=', bys, &
300                                      'byn=', byn, 'nx=', nx, 'ny=', ny
301             ENDIF
302             CALL local_stop
303          ENDIF
304
305!
306!--       Set the individual index arrays for all velocity components and
307!--       scalars, taking into account the staggered grid. The horizontal
308!--       wind component normal to a wall defines the position of the wall, and
309!--       in the respective direction the building is as long as specified in
310!--       building_length_?, but in the other horizontal direction (for w and s
311!--       in both horizontal directions) the building appears shortened by one
312!--       grid length due to the staggered grid.
313          nzb_local = 0
314          nzb_local(bys:byn-1,bxl:bxr-1) = bh
315
316       CASE ( 'read_from_file' )
317!
318!--       Arbitrary irregular topography data in PALM format (exactly matching
319!--       the grid size and total domain size)
320          OPEN( 90, FILE='TOPOGRAPHY_DATA', STATUS='OLD', FORM='FORMATTED',  &
321               ERR=10 )
322          DO  j = ny, 0, -1
323             READ( 90, *, ERR=11, END=11 )  ( topo_height(j,i), i = 0, nx )
324          ENDDO
325!
326!--       Calculate the index height of the topography
327          DO  i = 0, nx
328             DO  j = 0, ny
329                nzb_local(j,i) = NINT( topo_height(j,i) / dz )
330             ENDDO
331          ENDDO
332          nzb_local(-1,0:nx)   = nzb_local(ny,0:nx)
333          nzb_local(ny+1,0:nx) = nzb_local(0,0:nx)
334          nzb_local(:,-1)      = nzb_local(:,nx)
335          nzb_local(:,nx+1)    = nzb_local(:,0)
336     
337          GOTO 12
338
339 10       IF ( myid == 0 )  THEN
340             PRINT*, '+++ init_grid: file TOPOGRAPHY_DATA does not exist'
341          ENDIF
342          CALL local_stop
343
344 11       IF ( myid == 0 )  THEN
345             PRINT*, '+++ init_grid: errors in file TOPOGRAPHY_DATA'
346          ENDIF
347          CALL local_stop
348
349 12       CLOSE( 90 )
350
351       CASE DEFAULT
352!
353!--       The DEFAULT case is reached either if the parameter topography
354!--       contains a wrong character string or if the user has coded a special
355!--       case in the user interface. There, the subroutine user_init_grid
356!--       checks which of these two conditions applies.
357          CALL user_init_grid( nzb_local )
358
359    END SELECT
360
361!
362!-- Consistency checks and index array initialization are only required for
363!-- non-flat topography
364    IF ( TRIM( topography ) /= 'flat' )  THEN
365
366!
367!--    Consistency checks
368       IF ( MINVAL( nzb_local ) < 0  .OR.  MAXVAL( nzb_local ) > nz + 1 )  THEN
369          IF ( myid == 0 )  THEN
370             PRINT*, '+++ init_grid: nzb_local values are outside the', &
371                          'model domain'
372             PRINT*, '               MINVAL( nzb_local ) = ', MINVAL(nzb_local)
373             PRINT*, '               MAXVAL( nzb_local ) = ', MAXVAL(nzb_local)
374          ENDIF
375          CALL local_stop
376       ENDIF
377
378       IF ( bc_lr == 'cyclic' )  THEN
379          IF ( ANY( nzb_local(:,-1) /= nzb_local(:,nx)   )  .OR. &
380               ANY( nzb_local(:,0)  /= nzb_local(:,nx+1) ) )  THEN
381             IF ( myid == 0 )  THEN
382                PRINT*, '+++ init_grid: nzb_local does not fulfill cyclic', &
383                        '               boundary condition in x-direction'
384             ENDIF
385             CALL local_stop
386          ENDIF
387       ENDIF
388       IF ( bc_ns == 'cyclic' )  THEN
389          IF ( ANY( nzb_local(-1,:) /= nzb_local(ny,:)   )  .OR. &
390               ANY( nzb_local(0,:)  /= nzb_local(ny+1,:) ) )  THEN
391             IF ( myid == 0 )  THEN
392                PRINT*, '+++ init_grid: nzb_local does not fulfill cyclic', &
393                        '               boundary condition in y-direction'
394             ENDIF
395             CALL local_stop
396          ENDIF
397       ENDIF
398
399!
400!--    Initialize index arrays nzb_s_inner and nzb_w_inner
401       nzb_s_inner = nzb_local(nys-1:nyn+1,nxl-1:nxr+1)
402       nzb_w_inner = nzb_local(nys-1:nyn+1,nxl-1:nxr+1)
403
404!
405!--    Initialize remaining index arrays:
406!--    first pre-initialize them with nzb_s_inner...
407       nzb_u_inner = nzb_s_inner
408       nzb_u_outer = nzb_s_inner
409       nzb_v_inner = nzb_s_inner
410       nzb_v_outer = nzb_s_inner
411       nzb_w_outer = nzb_s_inner
412       nzb_s_outer = nzb_s_inner
413
414!
415!--    ...then extend pre-initialized arrays in their according directions
416!--    based on nzb_local using nzb_tmp as a temporary global index array
417
418!
419!--    nzb_s_outer:
420!--    extend nzb_local east-/westwards first, then north-/southwards
421       nzb_tmp = nzb_local
422       DO  j = -1, ny + 1
423          DO  i = 0, nx
424             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i), &
425                                 nzb_local(j,i+1) )
426          ENDDO
427       ENDDO
428       DO  i = nxl, nxr
429          DO  j = nys, nyn
430             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i), &
431                                     nzb_tmp(j+1,i) )
432          ENDDO
433!
434!--       non-cyclic boundary conditions (overwritten by call of
435!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
436          IF ( nys == 0 )  THEN
437             j = -1
438             nzb_s_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
439          ENDIF
440          IF ( nys == ny )  THEN
441             j = ny + 1
442             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
443          ENDIF
444       ENDDO
445!
446!--    nzb_w_outer:
447!--    identical to nzb_s_outer
448       nzb_w_outer = nzb_s_outer
449
450!
451!--    nzb_u_inner:
452!--    extend nzb_local rightwards only
453       nzb_tmp = nzb_local
454       DO  j = -1, ny + 1
455          DO  i = 0, nx + 1
456             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i) )
457          ENDDO
458       ENDDO
459       nzb_u_inner = nzb_tmp(nys-1:nyn+1,nxl-1:nxr+1)
460
461!
462!--    nzb_u_outer:
463!--    extend current nzb_tmp (nzb_u_inner) north-/southwards
464       DO  i = nxl, nxr
465          DO  j = nys, nyn
466             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i), &
467                                     nzb_tmp(j+1,i) )
468          ENDDO
469!
470!--       non-cyclic boundary conditions (overwritten by call of
471!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
472          IF ( nys == 0 )  THEN
473             j = -1
474             nzb_u_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
475          ENDIF
476          IF ( nys == ny )  THEN
477             j = ny + 1
478             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
479          ENDIF
480       ENDDO
481
482!
483!--    nzb_v_inner:
484!--    extend nzb_local northwards only
485       nzb_tmp = nzb_local
486       DO  i = -1, nx + 1
487          DO  j = 0, ny + 1
488             nzb_tmp(j,i) = MAX( nzb_local(j-1,i), nzb_local(j,i) )
489          ENDDO
490       ENDDO
491       nzb_v_inner = nzb_tmp(nys-1:nyn+1,nxl-1:nxr+1)
492
493!
494!--    nzb_v_outer:
495!--    extend current nzb_tmp (nzb_v_inner) right-/leftwards
496       DO  j = nys, nyn
497          DO  i = nxl, nxr
498             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i), &
499                                     nzb_tmp(j,i+1) )
500          ENDDO
501!
502!--       non-cyclic boundary conditions (overwritten by call of
503!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
504          IF ( nxl == 0 )  THEN
505             i = -1
506             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i+1), nzb_tmp(j,i) )
507          ENDIF
508          IF ( nxr == nx )  THEN
509             i = nx + 1
510             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i) )
511          ENDIF
512       ENDDO
513
514!
515!--    Exchange of lateral boundary values (parallel computers) and cyclic
516!--    boundary conditions, if applicable.
517!--    Since nzb_s_inner and nzb_w_inner are derived directly from nzb_local
518!--    they do not require exchange and are not included here.
519       CALL exchange_horiz_2d_int( nzb_u_inner )
520       CALL exchange_horiz_2d_int( nzb_u_outer )
521       CALL exchange_horiz_2d_int( nzb_v_inner )
522       CALL exchange_horiz_2d_int( nzb_v_outer )
523       CALL exchange_horiz_2d_int( nzb_w_outer )
524       CALL exchange_horiz_2d_int( nzb_s_outer )
525
526    ENDIF
527
528!
529!-- Preliminary: to be removed after completion of the topography code!
530!-- Set the former default k index arrays nzb_2d
531    nzb_2d      = nzb
532
533!
534!-- Set the individual index arrays which define the k index from which on
535!-- the usual finite difference form (which does not use surface fluxes) is
536!-- applied
537    IF ( prandtl_layer  .OR.  use_surface_fluxes )  THEN
538       nzb_diff_u         = nzb_u_inner + 2
539       nzb_diff_v         = nzb_v_inner + 2
540       nzb_diff_s_inner   = nzb_s_inner + 2
541       nzb_diff_s_outer   = nzb_s_outer + 2
542    ELSE
543       nzb_diff_u         = nzb_u_inner + 1
544       nzb_diff_v         = nzb_v_inner + 1
545       nzb_diff_s_inner   = nzb_s_inner + 1
546       nzb_diff_s_outer   = nzb_s_outer + 1
547    ENDIF
548
549!
550!-- Calculation of wall switches and factors required by diffusion_u/v.f90 and
551!-- for limitation of near-wall mixing length l_wall further below
552    corner_nl = 0
553    corner_nr = 0
554    corner_sl = 0
555    corner_sr = 0
556    wall_l    = 0
557    wall_n    = 0
558    wall_r    = 0
559    wall_s    = 0
560
561    DO  i = nxl, nxr
562       DO  j = nys, nyn
563!
564!--       u-component
565          IF ( nzb_u_outer(j,i) > nzb_u_outer(j+1,i) )  THEN
566             wall_u(j,i) = 1.0   ! north wall (location of adjacent fluid)
567             fym(j,i)    = 0.0
568             fyp(j,i)    = 1.0
569          ELSEIF ( nzb_u_outer(j,i) > nzb_u_outer(j-1,i) )  THEN
570             wall_u(j,i) = 1.0   ! south wall (location of adjacent fluid)
571             fym(j,i)    = 1.0
572             fyp(j,i)    = 0.0
573          ENDIF
574!
575!--       v-component
576          IF ( nzb_v_outer(j,i) > nzb_v_outer(j,i+1) )  THEN
577             wall_v(j,i) = 1.0   ! rigth wall (location of adjacent fluid)
578             fxm(j,i)    = 0.0
579             fxp(j,i)    = 1.0
580          ELSEIF ( nzb_v_outer(j,i) > nzb_v_outer(j,i-1) )  THEN
581             wall_v(j,i) = 1.0   ! left wall (location of adjacent fluid)
582             fxm(j,i)    = 1.0
583             fxp(j,i)    = 0.0
584          ENDIF
585!
586!--       w-component, also used for scalars, separate arrays for shear
587!--       production of tke
588          IF ( nzb_w_outer(j,i) > nzb_w_outer(j+1,i) )  THEN
589             wall_e_y(j,i) =  1.0   ! north wall (location of adjacent fluid)
590             wall_w_y(j,i) =  1.0
591             fwym(j,i)     =  0.0
592             fwyp(j,i)     =  1.0
593          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j-1,i) )  THEN
594             wall_e_y(j,i) = -1.0   ! south wall (location of adjacent fluid)
595             wall_w_y(j,i) =  1.0
596             fwym(j,i)     =  1.0
597             fwyp(j,i)     =  0.0
598          ENDIF
599          IF ( nzb_w_outer(j,i) > nzb_w_outer(j,i+1) )  THEN
600             wall_e_x(j,i) =  1.0   ! right wall (location of adjacent fluid)
601             wall_w_x(j,i) =  1.0
602             fwxm(j,i)     =  0.0
603             fwxp(j,i)     =  1.0
604          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j,i-1) )  THEN
605             wall_e_x(j,i) = -1.0   ! left wall (location of adjacent fluid)
606             wall_w_x(j,i) =  1.0
607             fwxm(j,i)     =  1.0
608             fwxp(j,i)     =  0.0
609          ENDIF
610!
611!--       Wall and corner locations inside buildings for limitation of
612!--       near-wall mixing length l_wall
613          IF ( nzb_s_inner(j,i) > nzb_s_inner(j+1,i) )  THEN
614
615             wall_n(j,i) = nzb_s_inner(j+1,i) + 1            ! North wall
616
617             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
618                corner_nl(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northleft corner
619                                      nzb_s_inner(j,i-1) ) + 1
620             ENDIF
621
622             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
623                corner_nr(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northright corner
624                                      nzb_s_inner(j,i+1) ) + 1
625             ENDIF
626
627          ENDIF
628
629          IF ( nzb_s_inner(j,i) > nzb_s_inner(j-1,i) )  THEN
630
631             wall_s(j,i) = nzb_s_inner(j-1,i) + 1            ! South wall
632             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
633                corner_sl(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southleft corner
634                                      nzb_s_inner(j,i-1) ) + 1
635             ENDIF
636
637             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
638                corner_sr(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southright corner
639                                      nzb_s_inner(j,i+1) ) + 1
640             ENDIF
641
642          ENDIF
643
644          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
645             wall_l(j,i) = nzb_s_inner(j,i-1) + 1            ! Left wall
646          ENDIF
647
648          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
649             wall_r(j,i) = nzb_s_inner(j,i+1) + 1            ! Right wall
650          ENDIF
651
652       ENDDO
653    ENDDO
654
655!
656!-- In case of topography: limit near-wall mixing length l_wall further:
657!-- Go through all points of the subdomain one by one and look for the closest
658!-- surface
659    IF ( TRIM(topography) /= 'flat' )  THEN
660       DO  i = nxl, nxr
661          DO  j = nys, nyn
662
663             nzb_si = nzb_s_inner(j,i)
664             vi     = vertical_influence(nzb_si)
665
666             IF ( wall_n(j,i) > 0 )  THEN
667!
668!--             North wall (y distance)
669                DO  k = wall_n(j,i), nzb_si
670                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i), 0.5 * dy )
671                ENDDO
672!
673!--             Above North wall (yz distance)
674                DO  k = nzb_si + 1, nzb_si + vi
675                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i),     &
676                                          SQRT( 0.25 * dy**2 + &
677                                          ( zu(k) - zw(nzb_si) )**2 ) )
678                ENDDO
679!
680!--             Northleft corner (xy distance)
681                IF ( corner_nl(j,i) > 0 )  THEN
682                   DO  k = corner_nl(j,i), nzb_si
683                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1), &
684                                               0.5 * SQRT( dx**2 + dy**2 ) )
685                   ENDDO
686!
687!--                Above Northleft corner (xyz distance)
688                   DO  k = nzb_si + 1, nzb_si + vi
689                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1),             &
690                                               SQRT( 0.25 * (dx**2 + dy**2) + &
691                                               ( zu(k) - zw(nzb_si) )**2 ) )
692                   ENDDO
693                ENDIF
694!
695!--             Northright corner (xy distance)
696                IF ( corner_nr(j,i) > 0 )  THEN
697                   DO  k = corner_nr(j,i), nzb_si
698                       l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1), &
699                                                0.5 * SQRT( dx**2 + dy**2 ) )
700                   ENDDO
701!
702!--                Above northright corner (xyz distance)
703                   DO  k = nzb_si + 1, nzb_si + vi
704                      l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1), &
705                                               SQRT( 0.25 * (dx**2 + dy**2) + &
706                                               ( zu(k) - zw(nzb_si) )**2 ) )
707                   ENDDO
708                ENDIF
709             ENDIF
710
711             IF ( wall_s(j,i) > 0 )  THEN
712!
713!--             South wall (y distance)
714                DO  k = wall_s(j,i), nzb_si
715                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i), 0.5 * dy )
716                ENDDO
717!
718!--             Above south wall (yz distance)
719                DO  k = nzb_si + 1, &
720                        nzb_si + vi
721                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i),     &
722                                          SQRT( 0.25 * dy**2 + &
723                                          ( zu(k) - zw(nzb_si) )**2 ) )
724                ENDDO
725!
726!--             Southleft corner (xy distance)
727                IF ( corner_sl(j,i) > 0 )  THEN
728                   DO  k = corner_sl(j,i), nzb_si
729                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1), &
730                                               0.5 * SQRT( dx**2 + dy**2 ) )
731                   ENDDO
732!
733!--                Above southleft corner (xyz distance)
734                   DO  k = nzb_si + 1, nzb_si + vi
735                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),             &
736                                               SQRT( 0.25 * (dx**2 + dy**2) + &
737                                               ( zu(k) - zw(nzb_si) )**2 ) )
738                   ENDDO
739                ENDIF
740!
741!--             Southright corner (xy distance)
742                IF ( corner_sr(j,i) > 0 )  THEN
743                   DO  k = corner_sr(j,i), nzb_si
744                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1), &
745                                               0.5 * SQRT( dx**2 + dy**2 ) )
746                   ENDDO
747!
748!--                Above southright corner (xyz distance)
749                   DO  k = nzb_si + 1, nzb_si + vi
750                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),             &
751                                               SQRT( 0.25 * (dx**2 + dy**2) + &
752                                               ( zu(k) - zw(nzb_si) )**2 ) )
753                   ENDDO
754                ENDIF
755
756             ENDIF
757
758             IF ( wall_l(j,i) > 0 )  THEN
759!
760!--             Left wall (x distance)
761                DO  k = wall_l(j,i), nzb_si
762                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1), 0.5 * dx )
763                ENDDO
764!
765!--             Above left wall (xz distance)
766                DO  k = nzb_si + 1, nzb_si + vi
767                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1),     &
768                                          SQRT( 0.25 * dx**2 + &
769                                          ( zu(k) - zw(nzb_si) )**2 ) )
770                ENDDO
771             ENDIF
772
773             IF ( wall_r(j,i) > 0 )  THEN
774!
775!--             Right wall (x distance)
776                DO  k = wall_r(j,i), nzb_si
777                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1), 0.5 * dx )
778                ENDDO
779!
780!--             Above right wall (xz distance)
781                DO  k = nzb_si + 1, nzb_si + vi
782                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1),     &
783                                          SQRT( 0.25 * dx**2 + &
784                                          ( zu(k) - zw(nzb_si) )**2 ) )
785                ENDDO
786
787             ENDIF
788
789          ENDDO
790       ENDDO
791
792    ENDIF
793
794!
795!-- Multiplication with wall_adjustment_factor
796    l_wall = wall_adjustment_factor * l_wall
797
798!
799!-- Need to set lateral boundary conditions for l_wall
800    CALL exchange_horiz( l_wall, 0, 0 )
801
802    DEALLOCATE( corner_nl, corner_nr, corner_sl, corner_sr, nzb_local, &
803                nzb_tmp, vertical_influence, wall_l, wall_n, wall_r, wall_s )
804
805
806 END SUBROUTINE init_grid
Note: See TracBrowser for help on using the repository browser.