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

Last change on this file since 39 was 39, checked in by raasch, 17 years ago

comments prepared for 3.1c

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