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