Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_tools / da_smooth_anl.inc
blob6a2bdf59a4dbfb1fc954b1455f370be68d46e567
1 subroutine da_smooth_anl(slab,imx,jmx,kx,npass,icrsdot)
3    !-----------------------------------------------------------------------
4    ! Purpose: spatially smooth (usually slab) to remove high
5    ! frequency waves
6    !-----------------------------------------------------------------------
8    implicit none
9    
10    real,    intent(inout) :: SLAB(:,:,:)
11    integer, intent(in)    :: imx, jmx, kx
12    integer, intent(in)    :: npass
13    integer, intent(in)    :: icrsdot
14    
15    real, allocatable :: SLABNEW(:,:)
16    real              :: XNU(1:2)
17    integer           :: ie, je, k 
18    integer           :: loop, n, i, j
20    if (trace_use) call da_trace_entry("da_smooth_anl")
21    
22    allocate (slabnew(imx,jmx))
24    ie=imx-1-icrsdot
25    je=jmx-1-icrsdot
26    xnu(1)=0.50
27    xnu(2)=-0.52
28    do k=1,kx
29       do loop=1,npass*2
30          n=2-mod(loop,2)
32          ! first smooth in the imx direction
34          do i=2,ie
35             do j=2,je
36                slabnew(i,j)=slab(i,j,k)+xnu(n) * &
37                ((slab(i,j+1,k)+slab(i,j-1,k))*0.5-slab(i,j,k))
38             end do
39          end do
40          do i=2,ie
41             do j=2,je
42                slab(i,j,k)=slabnew(i,j)
43             end do
44          end do
46          ! now smooth in the jmx direction
48          do j=2,je
49             do i=2,ie
50                slabnew(i,j)=slab(i,j,k)+xnu(n) * &
51                ((slab(i+1,j,k)+slab(i-1,j,k))*0.5-slab(i,j,k))
52             end do
53          end do
55          do i=2,ie
56             do j=2,je
57                slab(i,j,k)=slabnew(i,j)
58             end do
59          end do
60       end do
61    end do
63    deallocate (slabnew)
65    if (trace_use) call da_trace_exit("da_smooth_anl")
67 end subroutine da_smooth_anl