da_transform_through_rf_adj.inc

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