da_transform_through_rf.inc

References to this file elsewhere.
1 subroutine da_transform_through_rf(mz, grid_box_area, rf_alpha, val,&
2                                     xp, field)
3 
4    !---------------------------------------------------------------------------
5    ! Purpose: Control variable transform through the recursive filter.  
6    !
7    ! Method: 1) Transform to nondimensional v_hat space.
8    !         2) Perform recursive filter in non-dimensional space (ds = 1).
9    !         3) Scale recursive filter output.
10    !         4) Transform filtered field to dimensional space.
11    !         5) Optionally scale by background error.
12    !---------------------------------------------------------------------------
13 
14    implicit none
15 
16    type (xpose_type), intent(inout) :: xp   ! Dimensions and xpose buffers.
17    integer, intent(in) :: mz                ! Vertical truncation.
18 
19    real, intent(in)    :: grid_box_area(ims:ime,jms:jme) ! From first guess. 
20    real, intent(in)    :: rf_alpha(mz)      ! RF scale parameter. 
21    real, intent(in)    :: val(jds:jde,mz)   ! Error standard deviation.
22    real, intent(inout) :: field(ims:ime,jms:jme,kms:kme) ! Field to be transformed. 
23       
24    integer             :: rf_passes_over_two ! rf_passes / 2
25    integer             :: i, j, m, n, pass   ! Loop counters.
26    real                :: p_x(ims:ime,jms:jme)! sqrt(Grid box area).
27    real                :: val_j(xp%jtsy:xp%jtey)
28 
29    !----------------------------------------------------------------------
30    ! [1.0]: Initialise:
31    !----------------------------------------------------------------------  
32 
33    if (trace_use) call da_trace_entry("da_transform_through_rf")
34 
35    rf_passes_over_two = rf_passes / 2
36    
37    ! [1.1] Define inner product (square root of grid box area):
38    p_x(its:ite,jts:jte) = sqrt(grid_box_area(its:ite,jts:jte))
39 
40    ! [1.2] Transform to nondimensional v_hat space:
41 
42    do m = 1, mz
43       xp % v1z(its:ite,jts:jte,m) = field(its:ite,jts:jte,m) / &
44                                     p_x(its:ite,jts:jte)
45    end do
46 
47    !-------------------------------------------------------------------------
48    ! [2.0]: Perform 1D recursive filter in x-direction:
49    !-------------------------------------------------------------------------
50 
51    ! [2.1] Apply (i',j',k -> i,j',k') (xp % v1z -> xp % v1x)
52    ! convert from vertical column to x-stripe
53 
54    call da_transpose_z2x (xp)
55 
56    ! [2.2] Apply 1D filter in x direction:
57 
58    n = xp%itex - xp%itsx + 1
59    do m = xp%ktsx, min(xp%ktex,mz)
60       do j = xp%jtsx, xp%jtex
61          do pass = 1, rf_passes_over_two
62             call da_recursive_filter_1d(pass, rf_alpha(m), &
63                         xp % v1x(xp%itsx:xp%itex,j,m), n)
64          end do
65       end do
66    end do
67 
68    !-------------------------------------------------------------------------
69    ! [3.0]: Perform 1D recursive filter in y-direction:
70    !-------------------------------------------------------------------------
71 
72    ! [3.1] Apply (i, j' ,k' -> i', j ,k') (xp % v1x -> xp % v1y)
73    ! convert from vertical column to y-stripe
74 
75    call da_transpose_x2y (xp)
76 
77    ! [3.2] Apply 1D filter in y direction:
78 
79    n = xp%jtey - xp%jtsy + 1
80    do m = xp%ktsy, min(xp%ktey,mz)
81       do i = xp%itsy, xp%itey
82          val_j(xp%jtsy:xp%jtey) = xp%v1y(i,xp%jtsy:xp%jtey,m)
83          do pass = 1, rf_passes_over_two
84             call da_recursive_filter_1d(pass, rf_alpha(m), &
85                                          val_j, n)
86          end do
87          xp % v1y(i,xp%jtsy:xp%jtey,m) = val_j(xp%jtsy:xp%jtey)
88       end do
89    end do
90    
91    !-------------------------------------------------------------------------
92    ! [4.0]: Perform 1D recursive filter in y-direction:
93    !-------------------------------------------------------------------------
94 
95    ! [4.1] Apply (i',j,k' -> i',j',k) (xp % v1y -> xp % v1z)
96    ! convert from y-stripe to vertical column.
97    
98    call da_transpose_y2z (xp)
99 
100    ! [4.2] Transform filtered field to dimensional space:
101 
102    do m = 1, mz
103       field(its:ite,jts:jte,m) = xp % v1z(its:ite,jts:jte,m) *  &
104                                  p_x(its:ite,jts:jte)
105    end do
106 
107    ! [4.3] Optionally scale by background error:
108 
109    ! be_s % val = Gridpoint standard deviation - only required for
110    ! vert_corr = 1 as scaling is performed in vertical transform
111    ! for vert_corr = 2:
112 
113    if (vert_corr == 1) then
114       do m = 1, mz
115       do i = its, ite
116          field(i,jts:jte,m) = field(i,jts:jte,m) * val(jts:jte,m)
117       end do
118       end do
119    end if
120 
121    if (trace_use) call da_trace_exit("da_transform_through_rf")
122 
123 end subroutine da_transform_through_rf
124 
125