da_transform_vtovv.inc

References to this file elsewhere.
1 subroutine da_transform_vtovv(cv_size, xb, be, cv, vv, xp)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
6 
7    implicit none
8 
9    integer, intent(in)              :: cv_size ! Size of cv array.
10    type (xb_type), intent(in)       :: xb   ! First guess structure.
11    type (be_type), intent(in)       :: be   ! Background error structure.
12    real, intent(in)                 :: cv(1:cv_size)   ! control variables.
13    type (vp_type), intent(inout)    :: vv   ! Grid point/EOF equivalent.
14    type (xpose_type), intent(inout) :: xp   ! Dimensions and xpose buffers. 
15 
16    integer                          :: mz   ! Vertical truncation.
17 
18    if (trace_use) call da_trace_entry("da_transform_vtovv")
19 
20    !-------------------------------------------------------------------------
21    ! [1.0] Fill vv arrays from 1-dimensional cv array.
22    !-------------------------------------------------------------------------
23 
24    call da_cv_to_vv(cv_size, cv, xp, &
25       (/ be%v1%mz, be%v2%mz, be%v3%mz, be%v4%mz, be%v5%mz, be%ne /), vv)
26 
27    !-------------------------------------------------------------------------
28    ! [2.0] Perform VToVV Transform:
29    !-------------------------------------------------------------------------
30 
31    ! [2.1] Transform 1st control variable:
32 
33    mz = be % v1 % mz
34    if (mz > 0) then
35       call da_transform_through_rf(mz, xb % grid_box_area, be % v1 % rf_alpha,&
36                                     be % v1 % val, xp, vv % v1)
37    end if
38 
39    ! [2.2] Transform 2nd control variable:
40 
41    mz = be % v2 % mz
42    if (mz > 0) then
43       call da_transform_through_rf(mz, xb % grid_box_area, be % v2 % rf_alpha,&
44                                     be % v2 % val, xp, vv % v2)
45    end if
46 
47    ! [2.3] Transform 3rd control variable
48 
49    mz = be % v3 % mz
50    if (mz > 0) then
51       call da_transform_through_rf(mz, xb % grid_box_area, be % v3 % rf_alpha,&
52                                     be % v3 % val, xp, vv % v3)
53    end if
54 
55    ! [2.4] Transform 4th control variable
56       
57    mz = be % v4 % mz
58    if (mz > 0) then
59       call da_transform_through_rf(mz, xb % grid_box_area, be % v4 % rf_alpha,&
60                                     be % v4 % val, xp, vv % v4)
61    end if
62 
63    ! [2.5] Transform 5th control variable
64 
65    mz = be % v5 % mz
66    if (mz > 0) then
67       call da_transform_through_rf(mz, xb % grid_box_area, be % v5 % rf_alpha,&
68                                     be % v5 % val, xp, vv % v5)
69    end if
70 
71    ! [2.6] Transform alpha control variable
72 
73    mz = be % alpha % mz
74    if (mz > 0) then
75       call da_transform_through_rf(mz, xb % grid_box_area, be % alpha % rf_alpha,&
76                                     be % alpha % val, xp, vv % alpha)
77    end if
78 
79    if (trace_use) call da_trace_exit("da_transform_vtovv")
80 
81 
82 end subroutine da_transform_vtovv
83 
84