da_filter_adj.inc
References to this file elsewhere.
1 subroutine da_filter_adj(grid, var)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 !-----------------------------------------------------------------------
6
7 implicit none
8
9 type(domain), intent(inout) :: grid
10 real, intent(inout) :: var(ims:ime,jms:jme,kms:kme)
11
12 integer :: i, j, k
13
14 real :: w(3)
15
16 data w/0.25,0.5,0.25/
17
18 if (trace_use) call da_trace_entry("da_filter_adj")
19
20 ! Copy var for transpose.
21
22 grid%xp%v1z(its:ite,jts:jte,kts:kte) = var(its:ite,jts:jte,kts:kte)
23
24 ! Apply (i',j',k -> i',j,k') transpose (v1z -> v1y).
25 call da_transpose_z2y (grid)
26
27 ! Perform y-direction filter:
28 do k=grid%xp%ktsy, grid%xp%ktey
29 do i=grid%xp%itey, grid%xp%itsy, -1
30 ! Backward
31 do j=jds+1, jde-1
32 grid%xp%v1y(i,j-1,k) = grid%xp%v1y(i,j-1,k) + w(1)*grid%xp%v1y(i,j,k)
33 grid%xp%v1y(i,j+1,k) = grid%xp%v1y(i,j+1,k) + w(3)*grid%xp%v1y(i,j,k)
34 grid%xp%v1y(i,j ,k) = w(2)*grid%xp%v1y(i,j,k)
35 end do
36
37 ! Forward
38 do j=jde-1,jds+1,-1
39 grid%xp%v1y(i,j-1,k) = grid%xp%v1y(i,j-1,k) + w(1)*grid%xp%v1y(i,j,k)
40 grid%xp%v1y(i,j+1,k) = grid%xp%v1y(i,j+1,k) + w(3)*grid%xp%v1y(i,j,k)
41 grid%xp%v1y(i,j ,k) = w(2)*grid%xp%v1y(i,j,k)
42 end do
43 end do
44 end do
45
46 ! Apply (i',j,k' -> i,j',k') transpose (v1y -> v1x).
47 call da_transpose_y2x (grid)
48
49 ! Perform x-direction filter:
50 do k=grid%xp%ktsx, grid%xp%ktex
51 do j=grid%xp%jtex, grid%xp%jtsx, -1
52 ! Backward
53 do i=ids+1, ide-1
54 grid%xp%v1x(i-1,j,k) = grid%xp%v1x(i-1,j,k) + w(1)*grid%xp%v1x(i,j,k)
55 grid%xp%v1x(i+1,j,k) = grid%xp%v1x(i+1,j,k) + w(3)*grid%xp%v1x(i,j,k)
56 grid%xp%v1x(i,j ,k) = w(2)*grid%xp%v1x(i,j,k)
57 end do
58
59 ! Forward
60 do i=ide-1,ids+1,-1
61 grid%xp%v1x(i-1,j,k) = grid%xp%v1x(i-1,j,k) + w(1)*grid%xp%v1x(i,j,k)
62 grid%xp%v1x(i+1,j,k) = grid%xp%v1x(i+1,j,k) + w(3)*grid%xp%v1x(i,j,k)
63 grid%xp%v1x(i,j ,k) = w(2)*grid%xp%v1x(i,j,k)
64 end do
65 end do
66 end do
67
68 ! Apply (i,j',k' -> i',j',k) transpose (v1x -> v1z).
69 call da_transpose_x2z (grid)
70
71 var(its:ite,jts:jte,kts:kte) = grid%xp%v1z(its:ite,jts:jte,kts:kte)
72
73 if (trace_use) call da_trace_exit("da_filter_adj")
74
75 end subroutine da_filter_adj
76
77