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