wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_vtox_transforms / da_transform_vtovv_adj.inc
blob225b90dd66bdf660a618865c10b83a758b6c0b74
1 subroutine da_transform_vtovv_adj(grid, cv_size, be, cv, vv)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
9    type(domain),  intent(inout) :: grid
10    integer,       intent(in)    :: cv_size ! Size of cv array.
11    type(be_type), intent(in)    :: be   ! Background error structure.
12    real,          intent(inout) :: cv(1:cv_size)   ! control variables.
13    type(vp_type), intent(inout) :: vv   ! Grid point/EOF control var.
15    integer :: n    ! Loop counter.
16    integer :: mz   ! Vertical truncation.
17    integer :: ne   ! Ensemble size.
19    if (trace_use) call da_trace_entry("da_transform_vtovv_adj")
21    !-------------------------------------------------------------------------
22    ! [2.0] Perform VToVV Transform:
23    !-------------------------------------------------------------------------
25    ! [2.1] Transform 1st control variable:
27    mz = be % v1 % mz
28    if (mz > 0) then
29       call da_transform_through_rf_adj(grid, mz, be % v1 % rf_alpha, be % v1 % val, vv % v1)
30    end if
32    ! [2.2] Transform 2nd control variable:
34    mz = be % v2 % mz
35    if (mz > 0) then
36       call da_transform_through_rf_adj(grid, mz, be % v2 % rf_alpha, be % v2 % val, vv % v2)
37    end if
39    ! [2.3] Transform 3rd control variable
41    mz = be % v3 % mz
42    if (mz > 0) then
43       call da_transform_through_rf_adj(grid, mz, be % v3 % rf_alpha, be % v3 % val, vv % v3)
44    end if
45    
46    ! [2.4] Transform 4th control variable
47       
48    mz = be % v4 % mz
49    if (mz > 0) then
50       call da_transform_through_rf_adj(grid, mz, be % v4 % rf_alpha, be % v4 % val, vv % v4)
51    end if
53    ! [2.5] Transform 5th control variable
55    mz = be % v5 % mz
56    if (mz > 0) then
57       call da_transform_through_rf_adj(grid, mz, be % v5 % rf_alpha, be % v5 % val, vv % v5)
58    end if
60    ! [2.6] Transform alpha control variable
62    ne = be % ne
63    if (ne > 0) then
64       mz = be % alpha % mz
65       do n = 1, ne
66          call da_transform_through_rf_adj(grid, mz, be % alpha % rf_alpha, be % alpha % val, vv % alpha(:,:,:,n))
67       end do
68    end if
70    !-------------------------------------------------------------------------
71    ! [1.0] Fill 1D cv array from 3-dimensional vv arrays.
72    !-------------------------------------------------------------------------
74    call da_vv_to_cv( vv, grid%xp, (/ be%v1%mz, be%v2%mz, be%v3%mz, be%v4%mz, be%v5%mz, be%alpha%mz, be%ne /), &
75                      cv_size, cv)
77    if (trace_use) call da_trace_exit("da_transform_vtovv_adj")
79 end subroutine da_transform_vtovv_adj