da_vertical_transform.inc

References to this file elsewhere.
1 subroutine da_vertical_transform(string, be, vertical_wgt, vv, vp)
2 
3    !---------------------------------------------------------------------
4    ! Purpose: TBD
5    !---------------------------------------------------------------------
6 
7    implicit none   
8 
9    character(len=*), intent(in)    :: string      ! Character operation
10    type (be_type),   intent(in)    :: be          ! Background error structure.
11    real,             intent(in)    :: vertical_wgt(ims:ime,jms:jme,kms:kme) ! Weighting.
12    type (vp_type),   intent(inout) :: vv          ! CV in gridpt/EOF space.
13    type (vp_type),   intent(inout) :: vp          ! CV in gridpt/level space.
14 
15    integer                         :: j, m        ! Loop counters.
16    real                            :: alpha_stddev_inv ! 1/ sigma_alpha
17 
18    if (trace_use) call da_trace_entry("da_vertical_transform")
19 
20    select case(string)
21       
22    case ('u');
23       
24       !-------------------------------------------------------------------
25       ! [1.0] Perform vp(i,j,k) = E L^{1/2} vv(i,j,m) transform:
26       !------------------------------------------------------------------- 
27 
28       if (be % v1 % mz > 0) then
29          call da_transform_vvtovp (be % v1 % evec, be % v1 % val, vertical_wgt, &
30             vv % v1, vp % v1, be % v1 % mz, kte)
31       else
32          vp % v1(its:ite,jts:jte,kts:kte) = 0.0
33       end if
34 
35       if (be % v2 % mz > 0) then
36          call da_transform_vvtovp (be % v2 % evec, be % v2 % val, vertical_wgt, &
37             vv % v2, vp % v2, be % v2 % mz, kte)
38       else
39          vp % v2(its:ite,jts:jte,kts:kte) = 0.0
40       end if
41 
42       if (be % v3 % mz > 0) then
43          call da_transform_vvtovp (be % v3 % evec, be % v3 % val, vertical_wgt, &
44             vv % v3, vp % v3, be % v3 % mz, kte)
45       else
46          vp % v3(its:ite,jts:jte,kts:kte) = 0.0
47       end if
48 
49       if (be % v4 % mz > 0) then
50          call da_transform_vvtovp (be % v4 % evec, be % v4 % val, vertical_wgt, &
51             vv % v4, vp % v4, be % v4 % mz, kte)
52       else
53          vp % v4(its:ite,jts:jte,kts:kte) = 0.0
54       end if
55 
56       if (be % v5 % mz > 0) then
57          if (global) then
58             vp % v5(its:ite,jts:jte,1) = vv % v5(its:ite,jts:jte,1)
59          else 
60             call da_transform_vvtovp (be % v5 % evec, be % v5 % val, vertical_wgt, & 
61                vv % v5, vp % v5, be % v5 % mz, kts)    
62          end if
63       else
64          vp % v5(its:ite,jts:jte,kts:kts) = 0.0
65       end if
66 
67       if (be % ne > 0) then
68          do m = 1, be % ne
69             do j = jts, jte
70                vp % alpha(its:ite,j,m) = vv % alpha(its:ite,j,m) * be % alpha % val(j,m)
71             end do
72          end do
73       end if
74 
75    case ('u_inv');
76      
77       !------------------------------------------------------------------- 
78       ! [2.0] Perform vv(i,j,m) = L^{-1/2} E^T vp(i,j,k) transform:
79       !------------------------------------------------------------------- 
80 
81       if (be % v1 % mz > 0) then
82          call da_transform_vptovv (be % v1 % evec, be % v1 % val, vertical_wgt, &
83             vp % v1, vv % v1, be % v1 % mz, kds,kde, ims,ime, jms,jme, kms,kme, &
84             its,ite, jts,jte, kts,kte)
85       end if
86 
87       if (be % v2 % mz > 0) then
88          call da_transform_vptovv (be % v2 % evec, be % v2 % val, vertical_wgt, &
89             vp % v2, vv % v2, be % v2 % mz, kds,kde, ims,ime, jms,jme, kms,kme, &
90             its,ite, jts,jte, kts,kte)
91       end if
92 
93       if (be % v3 % mz > 0) then
94          call da_transform_vptovv (be % v3 % evec, be % v3 % val, vertical_wgt, &
95             vp % v3, vv % v3, be % v3 % mz, kds,kde, ims,ime, jms,jme, kms,kme, &
96             its,ite, jts,jte, kts,kte)
97       end if
98 
99       if (be % v4 % mz > 0) then
100          call da_transform_vptovv (be % v4 % evec, be % v4 % val, vertical_wgt, &
101             vp % v4, vv % v4, be % v4 % mz, kds,kde, ims,ime, jms,jme, kms,kme, &
102             its,ite, jts,jte, kts,kte)
103       end if
104 
105       if (be % v5 % mz > 0) then
106          if (global) then
107             vv % v5(its:ite,jts:jte,1) = vp % v5(its:ite,jts:jte,1)
108          else
109             call da_transform_vptovv (be % v5 % evec, be % v5 % val, vertical_wgt, &
110                vp % v5, vv % v5, be % v5 % mz, kds,kde, ims,ime, jms,jme, kms,kme, &
111                its,ite, jts,jte, kts,kte)
112          end if
113       end if
114 
115       if (be % ne > 0) then
116          do m = 1, be % ne
117             do j = jts, jte
118                alpha_stddev_inv = 1.0 / be % alpha % val(j,m)
119                vv % alpha(its:ite,j,m) = vp % alpha(its:ite,j,m) * alpha_stddev_inv
120             end do
121          end do
122       end if
123 
124    case ('u_adj');
125     
126       !------------------------------------------------------------------- 
127       ! [3.0] Perform vv_adj = U_{v}^{T} vp_adj transform:
128       !------------------------------------------------------------------- 
129 
130       if (be % v1 % mz > 0) then
131          call da_transform_vvtovp_adj (be % v1 % evec, be % v1 % val, vertical_wgt, &
132             vp % v1, vv % v1, be % v1 % mz, kte)
133       end if
134 
135       if (be % v2 % mz > 0) then
136          call da_transform_vvtovp_adj (be % v2 % evec, be % v2 % val, vertical_wgt, &
137             vp % v2, vv % v2, be % v2 % mz, kte)
138       end if
139 
140       if (be % v3 % mz > 0) then
141          call da_transform_vvtovp_adj (be % v3 % evec, be % v3 % val, vertical_wgt, &
142             vp % v3, vv % v3, be % v3 % mz, kte)
143       end if
144 
145       if (be % v4 % mz > 0) then
146          call da_transform_vvtovp_adj (be % v4 % evec, be % v4 % val, vertical_wgt, &
147             vp % v4, vv % v4, be % v4 % mz, kte)
148       end if
149 
150       if (be % v5 % mz > 0) then
151          if (global) then
152             vv % v5(its:ite,jts:jte,1) = vp % v5(its:ite,jts:jte,1)
153          else
154             call da_transform_vvtovp_adj (be % v5 % evec, be % v5 % val, vertical_wgt, &
155                vp % v5, vv % v5, be % v5 % mz, kts)
156          end if
157       end if
158 
159       if (be % ne > 0) then
160          do m = 1, be % ne
161             do j = jts, jte
162                vv % alpha(its:ite,j,m) = vp % alpha(its:ite,j,m) * be % alpha % val(j,m)
163             end do
164          end do
165       end if
166 
167    case default;
168    
169       call da_error(__FILE__,__LINE__, &
170          (/"Invalid da_vertical_transform option "//trim(string)/))
171 
172    end select
173 
174    if (trace_use) call da_trace_exit("da_vertical_transform")
175 
176 end subroutine da_vertical_transform
177 
178