da_check_vtox_adjoint.inc

References to this file elsewhere.
1 subroutine da_check_vtox_adjoint(grid, cv_size, xbx, be, ep, cv1, vv, vp)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Test V to X routine and adjoint for compatibility.
5    !
6    ! Method:  Standard adjoint test: < x, x > = < v_adj, v >.
7    !---------------------------------------------------------------------------
8 
9    implicit none
10 
11    type(domain), intent(inout)               :: grid
12 
13    integer, intent(in)             :: cv_size
14    type (xbx_type),intent(in)      :: xbx    ! For header & non-grid arrays.
15    type (be_type), intent(in)      :: be     ! background error structure.
16    type (ep_type), intent(in)      :: ep     ! Ensemble perturbation structure.
17    real, intent(in)                :: cv1(1:cv_size) ! control variable (local).
18    type (vp_type), intent(inout)   :: vv     ! Grdipt/EOF CV.
19    type (vp_type), intent(inout)   :: vp     ! Grdipt/level CV.
20 
21    real                   :: cv2(1:cv_size)    ! control variable (local).
22    real                   :: adj_par_lhs ! < x, x >
23    real                   :: adj_par_rhs ! < x, x >
24    real                   :: adj_sum_lhs ! < x, x >
25    real                   :: adj_sum_rhs ! < v_adj, v >
26 
27    if (trace_use) call da_trace_entry("da_check_vtox_adjoint")
28 
29    write(unit=stdout, fmt='(/a/)') &
30       'da_check_vtox_adjoint: Adjoint Test Results:'
31 
32    !----------------------------------------------------------------------
33    ! [1.0] Initialise:
34    !----------------------------------------------------------------------
35 
36    cv2(:) = 0.0
37       
38    !----------------------------------------------------------------------
39    ! [2.0] Perform x = U v transform:
40    !----------------------------------------------------------------------
41 
42    call da_zero_x (grid%xa)
43 
44    call da_transform_vtox(grid,cv_size, xbx, be, ep, cv1, vv, vp)
45 
46    !----------------------------------------------------------------------
47    ! [3.0] Calculate LHS of adjoint test equation: 
48    !----------------------------------------------------------------------
49 
50    adj_par_lhs = sum(grid%xa % u(its:ite, jts:jte, kts:kte)**2) / typical_u_rms**2 &
51                + sum(grid%xa % v(its:ite, jts:jte, kts:kte)**2) / typical_v_rms**2 &     
52                + sum(grid%xa % p(its:ite, jts:jte, kts:kte)**2) / typical_p_rms**2 &     
53                + sum(grid%xa % t(its:ite, jts:jte, kts:kte)**2) / typical_t_rms**2 &     
54                + sum(grid%xa % q(its:ite, jts:jte, kts:kte)**2) / typical_q_rms**2 &     
55                + sum(grid%xa % rho(its:ite,jts:jte,kts:kte)**2)/ typical_rho_rms**2 & 
56                + sum(grid%xa % psfc(its:ite, jts:jte)**2) / typical_p_rms**2             
57 
58    if (cv_options_hum == cv_options_hum_relative_humidity) then
59       adj_par_lhs = adj_par_lhs &
60               + sum(grid%xa % rh(its:ite, jts:jte, kts:kte)**2) / typical_rh_rms**2
61    end if
62 
63    if (use_radarobs) then
64       adj_par_lhs = adj_par_lhs &
65          + sum(grid%xa % wh (its:ite, jts:jte, kts:kte)**2)/typical_w_rms**2
66    else
67       adj_par_lhs = adj_par_lhs &
68          + sum(grid%xa % w  (its:ite, jts:jte, kts:kte)**2)/typical_w_rms**2
69    end if
70 
71    !-------------------------------------------------------------------------
72    ! [4.0] Rescale input to adjoint routine:
73    !-------------------------------------------------------------------------
74 
75    grid%xa % u(:,:,:) = grid%xa % u(:,:,:) / typical_u_rms**2
76    grid%xa % v(:,:,:) = grid%xa % v(:,:,:) / typical_v_rms**2
77    grid%xa % p(:,:,:) = grid%xa % p(:,:,:) / typical_p_rms**2
78    grid%xa % t(:,:,:) = grid%xa % t(:,:,:) / typical_t_rms**2
79    grid%xa % q(:,:,:) = grid%xa % q(:,:,:) / typical_q_rms**2
80    grid%xa % rho(:,:,:) = grid%xa % rho(:,:,:) / typical_rho_rms**2
81 
82    grid%xa % psfc(:,:) = grid%xa % psfc(:,:) / typical_p_rms**2
83 
84    if (cv_options_hum == cv_options_hum_relative_humidity) then
85       grid%xa % rh(:,:,:) = grid%xa % rh(:,:,:) / typical_rh_rms**2
86    end if
87 
88    if (use_radarobs) then
89       grid%xa %wh(:,:,:) = grid%xa %wh(:,:,:) / typical_w_rms**2
90       grid%xa % w(:,:,:) = 0.0
91    else
92       grid%xa %w (:,:,:) = grid%xa %w (:,:,:) / typical_w_rms**2
93    end if
94 
95    !-------------------------------------------------------------------------
96    ! [5.0] Perform adjoint operation:
97    !-------------------------------------------------------------------------
98 
99    call da_transform_vtox_adj(grid, cv_size, xbx, be, ep, vp, vv, cv2)
100 
101    !-------------------------------------------------------------------------
102    ! [6.0] Calculate RHS of adjoint test equation:
103    !-------------------------------------------------------------------------
104 
105    adj_par_rhs = sum(cv1(1:cv_size) * cv2(1:cv_size))
106 
107    !-------------------------------------------------------------------------
108    ! [7.0] Print output:
109    !-------------------------------------------------------------------------
110 
111    adj_sum_lhs = wrf_dm_sum_real(adj_par_lhs)
112 
113    if (global) then
114       adj_sum_rhs = adj_par_rhs
115    else
116       adj_sum_rhs = wrf_dm_sum_real(adj_par_rhs)
117    end if
118    write (unit=stdout,fmt='(A,2F15.2)') &
119       'TEST_COVERAGE_da_check_vtox_adjoint:  adj_sum_lhs,adj_sum_rhs = ', &
120       adj_sum_lhs,adj_sum_rhs
121 
122    if (rootproc) then
123       write(unit=stdout, fmt='(/)')
124       write(unit=stdout, fmt='(a,1pe22.14)') &
125          'Whole  Domain: < x, x >     = ', adj_sum_lhs, &
126          'Whole  Domain: < v_adj, v > = ', adj_sum_rhs
127    end if
128 
129    write(unit=stdout, fmt='(/a/)') 'da_check_vtox_adjoint: Finished'
130 
131    if (trace_use) call da_trace_exit("da_check_vtox_adjoint")
132 
133 end subroutine da_check_vtox_adjoint
134 
135