da_test_vtoy_transform.inc

References to this file elsewhere.
1 subroutine da_test_vtoy_transform(grid, config_flags, vp, vv, &
2                                   xb, xbx, xp, be, xa, iv, y)
3 
4    !---------------------------------------------------------------------------
5    !  Purpose: Perform V To Y Adjoint transform test                             
6    !
7    !  Method:  Perform adjoint test on complete transform: <y,y> = <v_adj,v>.
8    !---------------------------------------------------------------------------
9 
10    implicit none
11 
12    type(grid_config_rec_type), intent(inout) :: config_flags
13    type(domain) ,   target              :: grid
14    type (vp_type),    intent(inout)     :: vv    ! Grdipt/EOF CV.
15    type (vp_type),    intent(inout)     :: vp    ! Grdipt/level CV.
16    type (xb_type),    intent(in)        :: xb    ! first guess (local).
17    type (xbx_type),   intent(in)        :: xbx   ! Header & non-gridded vars.
18    type (xpose_type), intent(inout)     :: xp    ! Dimensions and xpose buffers.
19    type (be_type),    intent(in)        :: be    ! background error structure.
20    type (x_type),     intent(out)       :: xa    ! analysis increments (local).
21    type (ob_type),    intent(inout)     :: iv    ! ob. increment vector.
22    type (y_type),     intent(inout)     :: y     ! y = h (xa)
23 
24 
25    real, dimension(1:cv_size)        :: cv          ! Test control variable.
26    real, dimension(1:cv_size)        :: cv_2
27 
28    real, dimension(ims:ime, jms:jme) :: field
29 
30    real                           :: adj_sum_lhs               ! < y, y >
31    real                           :: adj_rhs,adj_sum_rhs       ! < v, v_adj >
32    real                           :: partial_lhs   ! < y, y >
33    real                           :: pertile_lhs   ! < y, y >
34 
35    write(unit=*, fmt='(/a/a)') &
36         '       da_test_vtoy_transform:',&
37         '---------------------------------------'
38 
39 
40    call random_number(cv(:))
41    cv(:) = cv(:) - 0.5
42 
43    cv_2(1:cv_size) = cv(1:cv_size)
44 
45    call da_zero_x(xa)
46    call da_zero_vp_type(vp)
47    call da_zero_vp_type(vv)
48 
49    call da_transform_vtoy(be, cv, iv, vp, vv, xa, xb, xbx, xp, y, &
50                              grid, config_flags )
51 
52    !-------------------------------------------------------------------------
53    ! [3.0] Calculate LHS of adjoint test equation and
54    !       Rescale input to adjoint routine :
55    !-------------------------------------------------------------------------
56 
57    call da_get_y_lhs_value( iv, y, partial_lhs, pertile_lhs, adj_sum_lhs )
58 
59    cv = 0.0
60 
61    ! call da_zero_vp_type(vp)
62    ! call da_zero_vp_type(vv)
63    ! call da_zero_x(xa)      
64 
65    call da_transform_vtoy_adj(be, cv, iv, vp, vv, xa, xb, xbx, xp, y, &
66                                  grid, config_flags )
67 
68    adj_rhs = sum( cv(1:cv_size) * cv_2(1:cv_size) )
69    !-------------------------------------------------------------------------
70    ! Print output:
71    !-------------------------------------------------------------------------
72    ! JRB fix with wrf_dm_sum_real
73 
74 #ifdef DM_PARALLEL
75   if( global ) then
76      adj_sum_rhs = adj_rhs
77   else
78      call mpi_allreduce( adj_rhs, adj_sum_rhs, 1, true_mpi_real, mpi_sum, &
79                        comm, ierr)
80   end if
81 #else
82    adj_sum_rhs = adj_rhs
83    adj_sum_lhs = partial_lhs
84 #endif
85 
86 #ifdef DM_PARALLEL
87    if ( rootproc ) then
88       write(unit=*, fmt='(A,1pe22.14)') &
89          'Whole Domain  < y, y     > = ', adj_sum_lhs
90       write(unit=*, fmt='(A,1pe22.14)') &
91          'Whole Domain  < v, v_adj > = ', adj_sum_rhs
92    end if
93 #else
94    write(unit=*, fmt='(A,1pe22.14)') &
95       'Whole Domain  < y, y     > = ', adj_sum_lhs
96    write(unit=*, fmt='(A,1pe22.14)') &
97       'Whole Domain  < v, v_adj > = ', adj_sum_rhs
98 #endif
99 
100 
101 end subroutine da_test_vtoy_transform
102 
103