da_transform_through_rf_adj.inc

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