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

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

Initial repository layout and content

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