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 rf_passes_over_two = rf_passes / 2
34
35 ! [1.1] Define inner product (square root of grid box area):
36 p_x(its:ite,jts:jte) = sqrt(grid_box_area(its:ite,jts:jte))
37
38 ! [1.2] Transform to nondimensional v_hat space:
39
40 do m = 1, mz
41 xp % v1z(its:ite,jts:jte,m) = field(its:ite,jts:jte,m) / &
42 p_x(its:ite,jts:jte)
43 end do
44
45 !-------------------------------------------------------------------------
46 ! [2.0]: Perform 1D recursive filter in x-direction:
47 !-------------------------------------------------------------------------
48
49 ! [2.1] Apply (i',j',k -> i,j',k') (xp % v1z -> xp % v1x)
50 ! convert from vertical column to x-stripe
51
52 call da_transpose_z2x (xp)
53
54 ! [2.2] Apply 1D filter in x direction:
55
56 n = xp%itex - xp%itsx + 1
57 do m = xp%ktsx, min(xp%ktex,mz)
58 do j = xp%jtsx, xp%jtex
59 do pass = 1, rf_passes_over_two
60 call da_recursive_filter_1d(pass, rf_alpha(m), &
61 xp % v1x(xp%itsx:xp%itex,j,m), n)
62 end do
63 end do
64 end do
65
66 !-------------------------------------------------------------------------
67 ! [3.0]: Perform 1D recursive filter in y-direction:
68 !-------------------------------------------------------------------------
69
70 ! [3.1] Apply (i, j' ,k' -> i', j ,k') (xp % v1x -> xp % v1y)
71 ! convert from vertical column to y-stripe
72
73 call da_transpose_x2y (xp)
74
75 ! [3.2] Apply 1D filter in y direction:
76
77 n = xp%jtey - xp%jtsy + 1
78 do m = xp%ktsy, min(xp%ktey,mz)
79 do i = xp%itsy, xp%itey
80 val_j(xp%jtsy:xp%jtey) = xp%v1y(i,xp%jtsy:xp%jtey,m)
81 do pass = 1, rf_passes_over_two
82 call da_recursive_filter_1d(pass, rf_alpha(m), &
83 val_j, n)
84 end do
85 xp % v1y(i,xp%jtsy:xp%jtey,m) = val_j(xp%jtsy:xp%jtey)
86 end do
87 end do
88
89 !-------------------------------------------------------------------------
90 ! [4.0]: Perform 1D recursive filter in y-direction:
91 !-------------------------------------------------------------------------
92
93 ! [4.1] Apply (i',j,k' -> i',j',k) (xp % v1y -> xp % v1z)
94 ! convert from y-stripe to vertical column.
95
96 call da_transpose_y2z (xp)
97
98 ! [4.2] Transform filtered field to dimensional space:
99
100 do m = 1, mz
101 field(its:ite,jts:jte,m) = xp % v1z(its:ite,jts:jte,m) * &
102 p_x(its:ite,jts:jte)
103 end do
104
105 ! [4.3] Optionally scale by background error:
106
107 ! be_s % val = Gridpoint standard deviation - only required for
108 ! vert_corr = 1 as scaling is performed in vertical transform
109 ! for vert_corr = 2:
110
111 if (vert_corr == 1) then
112 do m = 1, mz
113 do i = its, ite
114 field(i,jts:jte,m) = field(i,jts:jte,m) * val(jts:jte,m)
115 end do
116 end do
117 end if
118
119 end subroutine da_transform_through_rf
120
121