da_check_vvtovp_adjoint.inc

References to this file elsewhere.
1 subroutine da_check_vvtovp_adjoint(ne, xb, be, vv, vp)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Test Vv to Vp routine and adjoint for compatibility.
5    !
6    ! Method:  Standard adjoint test: < Vp, Vp > = < Vv_adj, Vv >.
7    !---------------------------------------------------------------------------
8 
9    implicit none
10 
11    integer, intent(in)               :: ne    ! Ensemble size.
12    type (xb_type), intent(in)        :: xb    ! first guess (local).
13    type (be_type), intent(in)        :: be    ! background error structure.
14    type (vp_type), intent(inout)     :: vv    ! CV(i,j,m).
15    type (vp_type), intent(inout)     :: vp    ! CV(i,j,k)
16 
17    real                              :: adj_par_lhs ! < x, x >
18    real                              :: adj_par_rhs ! < v_adj, v >
19    real                              :: adj_sum_lhs ! < x, x >
20    real                              :: adj_sum_rhs ! < v_adj, v >
21 
22    real                              :: vv2_v1(ims:ime,jms:jme,kms:kme)
23    real                              :: vv2_v2(ims:ime,jms:jme,kms:kme)
24    real                              :: vv2_v3(ims:ime,jms:jme,kms:kme)
25    real                              :: vv2_v4(ims:ime,jms:jme,kms:kme)
26    real                              :: vv2_v5(ims:ime,jms:jme,kms:kme)
27    real                              :: vv2_alpha(ims:ime,jms:jme,1:ne)
28 
29    !----------------------------------------------------------------------
30    ! [1.0] Initialise:
31    !----------------------------------------------------------------------
32 
33    write(unit=*, fmt='(/a/)') 'da_check_vvtovp_adjoint: Test Results:'
34       
35    !----------------------------------------------------------------------
36    ! [2.0] Perform Vp = U_v Vv transform:
37    !----------------------------------------------------------------------
38 
39    call da_vertical_transform('u', be, &
40                                xb % vertical_inner_product, &
41                                vv, vp)
42 
43    !----------------------------------------------------------------------
44    ! [3.0] Calculate LHS of adjoint test equation:
45    !----------------------------------------------------------------------
46 
47    adj_par_lhs = sum(vp % v1(its:ite,jts:jte,kts:kte)**2) * inv_typ_vp1_sumsq &
48                + sum(vp % v2(its:ite,jts:jte,kts:kte)**2) * inv_typ_vp2_sumsq &
49                + sum(vp % v3(its:ite,jts:jte,kts:kte)**2) * inv_typ_vp3_sumsq &
50                + sum(vp % v4(its:ite,jts:jte,kts:kte)**2) * inv_typ_vp4_sumsq &
51                + sum(vp % v5(its:ite,jts:jte,kts:kte)**2) * inv_typ_vp5_sumsq
52 
53    if (be % ne > 0) then
54       adj_par_lhs = adj_par_lhs + &
55          sum(vp % alpha(its:ite,jts:jte,1:be%ne)**2) * inv_typ_vpalpha_sumsq
56    end if
57 
58    !----------------------------------------------------------------------
59    ! [4.0] Rescale input to adjoint routine:
60    !----------------------------------------------------------------------
61 
62    vp % v1(its:ite,jts:jte,kts:kte) = vp % v1(its:ite,jts:jte,kts:kte) * &
63       inv_typ_vp1_sumsq
64    vp % v2(its:ite,jts:jte,kts:kte) = vp % v2(its:ite,jts:jte,kts:kte) * &
65       inv_typ_vp2_sumsq
66    vp % v3(its:ite,jts:jte,kts:kte) = vp % v3(its:ite,jts:jte,kts:kte) * &
67       inv_typ_vp3_sumsq
68    vp % v4(its:ite,jts:jte,kts:kte) = vp % v4(its:ite,jts:jte,kts:kte) * &
69       inv_typ_vp4_sumsq
70    vp % v5(its:ite,jts:jte,kts:kte) = vp % v5(its:ite,jts:jte,kts:kte) * &
71       inv_typ_vp5_sumsq
72 
73    if (be % ne > 0) then
74       vp % alpha(its:ite,jts:jte,1:be%ne) = &
75          vp % alpha(its:ite,jts:jte,1:be%ne) * inv_typ_vpalpha_sumsq
76    end if
77 
78    !----------------------------------------------------------------------
79    ! [5.0] Perform adjoint operation:
80    !----------------------------------------------------------------------
81 
82    vv2_v1(its:ite,jts:jte,1:be%v1%mz) = vv % v1(its:ite,jts:jte,1:be%v1%mz)
83    vv2_v2(its:ite,jts:jte,1:be%v2%mz) = vv % v2(its:ite,jts:jte,1:be%v2%mz)
84    vv2_v3(its:ite,jts:jte,1:be%v3%mz) = vv % v3(its:ite,jts:jte,1:be%v3%mz)
85    vv2_v4(its:ite,jts:jte,1:be%v4%mz) = vv % v4(its:ite,jts:jte,1:be%v4%mz)
86    vv2_v5(its:ite,jts:jte,1:be%v5%mz) = vv % v5(its:ite,jts:jte,1:be%v5%mz)
87 
88    if (be % ne > 0) then
89       vv2_alpha(its:ite,jts:jte,1:be%ne) = vv % alpha(its:ite,jts:jte,1:be%ne)
90    end if      
91 
92    call da_vertical_transform('u_adj', be, &
93                                xb % vertical_inner_product, &
94                                vv, vp)
95 
96    !----------------------------------------------------------------------
97    ! [6.0] Calculate RHS of adjoint test equation:
98    !----------------------------------------------------------------------
99 
100    adj_par_rhs = 0.0
101    if (be % v1 % mz > 0) &
102       adj_par_rhs = sum(vv % v1(its:ite,jts:jte,1:be%v1%mz) * &
103          vv2_v1(its:ite,jts:jte,1:be%v1%mz)) + adj_par_rhs
104    if (be % v2 % mz > 0) &
105       adj_par_rhs = sum(vv % v2(its:ite,jts:jte,1:be%v2%mz) * &
106          vv2_v2(its:ite,jts:jte,1:be%v2%mz)) + adj_par_rhs
107    if (be % v3 % mz > 0) &
108       adj_par_rhs = sum(vv % v3(its:ite,jts:jte,1:be%v3%mz) * &
109          vv2_v3(its:ite,jts:jte,1:be%v3%mz)) + adj_par_rhs
110    if (be % v4 % mz > 0) &
111       adj_par_rhs = sum(vv % v4(its:ite,jts:jte,1:be%v4%mz) * &
112          vv2_v4(its:ite,jts:jte,1:be%v4%mz)) + adj_par_rhs
113    if (be % v5 % mz == 1) &
114       adj_par_rhs = sum(vv % v5(its:ite,jts:jte,1:be%v5%mz) * &
115          vv2_v5(its:ite,jts:jte,1:be%v5%mz)) + adj_par_rhs
116    if (be % ne > 0) &
117       adj_par_rhs = sum(vv % alpha(its:ite,jts:jte,1:be%ne) * &
118          vv2_alpha(its:ite,jts:jte,1:be%ne)) + adj_par_rhs
119 
120    !----------------------------------------------------------------------
121    ! [7.0] Print output:
122    !----------------------------------------------------------------------
123 
124    write(unit=stdout, fmt='(a,1pe22.14)') &
125         'Single domain < vp,     vp > = ', adj_par_lhs, &
126         'Single domain < Vv_adj, Vv > = ', adj_par_rhs
127 
128    adj_sum_lhs = wrf_dm_sum_real(adj_par_lhs)
129    adj_sum_rhs = wrf_dm_sum_real(adj_par_rhs)
130 
131    write(unit=stdout,fmt='(A,2F10.2)') &
132       'TEST_COVERAGE_da_check_vvtovp_adjoint:  adj_sum_lhs,adj_sum_rhs = ', &
133       adj_sum_lhs,adj_sum_rhs
134 
135    if (rootproc) then
136       write(unit=stdout, fmt='(/)')
137       write(unit=stdout, fmt='(a,1pe22.14)') &
138          'Whole  Domain: < Vp, Vp >     = ', adj_sum_lhs, &
139          'Whole  Domain: < Vv_adj, Vv > = ', adj_sum_rhs
140    end if
141       
142    vv % v1(its:ite,jts:jte,1:be%v1%mz) = vv2_v1(its:ite,jts:jte,1:be%v1%mz)
143    vv % v2(its:ite,jts:jte,1:be%v2%mz) = vv2_v2(its:ite,jts:jte,1:be%v2%mz)
144    vv % v3(its:ite,jts:jte,1:be%v3%mz) = vv2_v3(its:ite,jts:jte,1:be%v3%mz)
145    vv % v4(its:ite,jts:jte,1:be%v4%mz) = vv2_v4(its:ite,jts:jte,1:be%v4%mz)
146    vv % v5(its:ite,jts:jte,1:be%v5%mz) = vv2_v5(its:ite,jts:jte,1:be%v5%mz)
147 
148    if (be % ne > 0) then
149       vv % alpha(its:ite,jts:jte,1:be%ne) = vv2_alpha(its:ite,jts:jte,1:be%ne)
150    end if
151 
152    write(unit=stdout, fmt='(/a/)') 'da_check_vvtovp_adjoint: Test Finished.'
153 
154 end subroutine da_check_vvtovp_adjoint
155 
156