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

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

updating comments and rc-file

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