Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_bogus / da_transform_xtoy_bogus_adj.inc
blob777d87f71cfcfe131a23d2b4eac896c78d1b331a
1 subroutine da_transform_xtoy_bogus_adj(grid, iv, jo_grad_y, jo_grad_x)
3    !---------------------------------------------------------------------
4    ! Purpose: the adjoint of bogus observation operators.
5    !    Updated for Analysis on Arakawa-C grid
6    !    Author: Syed RH Rizvi,  MMM/ESSL/NCAR,  Date: 10/22/2008
7    !---------------------------------------------------------------------
9    implicit none
11    type (domain),     intent(in)     :: grid
12    type (iv_type),    intent(in)     :: iv          ! obs. inc vector (o-b).
13    type (y_type) ,    intent(inout)  :: jo_grad_y   ! grad_y(jo)
14    type (x_type) ,    intent(inout)  :: jo_grad_x   ! grad_x(jo)
16    integer :: n        ! Loop counter.
17    integer :: i, j, k  ! Index dimension.
18    integer :: num_levs
19    real    :: dx, dxm
20    real    :: dy, dym
22    real    :: model_t(kms:kme)
23    real    :: tt(kms:kme)
24    real    :: model_q(kms:kme)
25    real    :: qq(kms:kme)
26    real    :: model_p_c(kms:kme)
27    real    :: pp(kms:kme)
28    real    :: psfc_loc,terr_loc,ppsfc   
30    real, allocatable :: temp_u(:,:)
31    real, allocatable :: temp_v(:,:)
32    real, allocatable :: temp_t(:,:)
33    real, allocatable :: temp_q(:,:)           
35    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_bogus_adj")
37    allocate (temp_u(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
38    allocate (temp_v(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
39    allocate (temp_t(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
40    allocate (temp_q(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
42    do n=iv%info(bogus)%n1,iv%info(bogus)%n2
43       temp_u(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%u
44       temp_v(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%v
45       temp_t(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%t
46       temp_q(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%q
47    end do
49    ! [1.1] Adjoint feedback from Y to X for u and v:
51 #ifdef A2C
52    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(bogus), temp_u,'u')
53    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(bogus), temp_v,'v')
54 #else
55    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(bogus), temp_u)
56    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(bogus), temp_v)
57 #endif
58    call da_interp_lin_3d_adj (jo_grad_x%t, iv%info(bogus), temp_t)
59    call da_interp_lin_3d_adj (jo_grad_x%q, iv%info(bogus), temp_q)
60    deallocate (temp_u)
61    deallocate (temp_v)
62    deallocate (temp_t)
63    deallocate (temp_q)
65    do n= iv%info(bogus)%n1, iv%info(bogus)%n2
66       num_levs = iv%info(bogus)%levels(n)
67       if (num_levs < 1) cycle
69       ! [1.1] Get cross pt. horizontal interpolation weights:
71       i   = iv%info(bogus)%i(1,n)
72       dy  = iv%info(bogus)%dy(1,n)
73       dym = iv%info(bogus)%dym(1,n)
74       j   = iv%info(bogus)%j(1,n)
75       dx  = iv%info(bogus)%dx(1,n)
76       dxm = iv%info(bogus)%dxm(1,n)
78       ! [1.2] Compute the feedback from SLP to t and q:
80       ! 1.2.1 Background at the observation location:
82       do k = kts, kte
83          model_t(k) = dym*(dxm*grid%xb%t(i,j,k) + dx*grid%xb%t(i+1,j,k)) &
84             + dy *(dxm*grid%xb%t(i,j+1,k) + dx*grid%xb%t(i+1,j+1,k))
85          model_q(k) = dym*(dxm*grid%xb%q(i,j,k) + dx*grid%xb%q(i+1,j,k)) &
86             + dy *(dxm*grid%xb%q(i,j+1,k) + dx*grid%xb%q(i+1,j+1,k))
87          model_p_c(k) = dym*(dxm*grid%xb%p(i,j,k) + dx*grid%xb%p(i+1,j,k)) &
88             + dy *(dxm*grid%xb%p(i,j+1,k) + dx*grid%xb%p(i+1,j+1,k))
89       end do
91       terr_loc = dym*(dxm*grid%xb%terr(i,j)   + dx*grid%xb%terr(i+1,j)) &
92          + dy *(dxm*grid%xb%terr(i,j+1) + dx*grid%xb%terr(i+1,j+1))
93       psfc_loc = dym*(dxm*grid%xb%psfc(i,j)   + dx*grid%xb%psfc(i+1,j)) &
94          + dy *(dxm*grid%xb%psfc(i,j+1) + dx*grid%xb%psfc(i+1,j+1))
96       ! 1.2.2 Compute the feedback from SLP to p, t, q, and psfc 
97       !       at the observed location:
99       call da_tpq_to_slp_adj(model_t, model_q ,model_p_c, terr_loc, psfc_loc, &
100          tt, qq, pp, ppsfc, jo_grad_y%bogus(n)%slp)  
102       ! 1.2.3 feedback from the observed location to grid space:
104       ! 1.2.3.1 for Psfc
106       jo_grad_x % psfc(i,j)     = jo_grad_x % psfc(i,j) + dym*dxm*ppsfc
107       jo_grad_x % psfc(i+1,j)   = jo_grad_x % psfc(i+1,j)+ dym*dx*ppsfc
108       jo_grad_x % psfc(i,j+1)   = jo_grad_x % psfc(i,j+1)+ dy*dxm*ppsfc
109       jo_grad_x % psfc(i+1,j+1) = jo_grad_x % psfc(i+1,j+1)+dy*dx*ppsfc
111       ! 1.2.3.2 for t, q, p
113       do k = kts, kte
114          jo_grad_x % t(i,j,k)     = jo_grad_x % t(i,j,k)+dym*dxm*tt(k)
115          jo_grad_x % t(i+1,j,k)   = jo_grad_x % t(i+1,j,k)+ dym*dx*tt(k)
116          jo_grad_x % t(i,j+1,k)   = jo_grad_x % t(i,j+1,k)+ dy*dxm*tt(k)
117          jo_grad_x % t(i+1,j+1,k) = jo_grad_x % t(i+1,j+1,k)+dy*dx*tt(k)
118          jo_grad_x % q(i,j,k)     = jo_grad_x % q(i,j,k)+dym*dxm*qq(k)
119          jo_grad_x % q(i+1,j,k)   = jo_grad_x % q(i+1,j,k)+ dym*dx*qq(k)
120          jo_grad_x % q(i,j+1,k)   = jo_grad_x % q(i,j+1,k)+ dy*dxm*qq(k)
121          jo_grad_x % q(i+1,j+1,k) = jo_grad_x % q(i+1,j+1,k)+dy*dx*qq(k)
122          jo_grad_x % p(i,j,k)     = jo_grad_x % p(i,j,k)+dym*dxm*pp(k)
123          jo_grad_x % p(i+1,j,k)   = jo_grad_x % p(i+1,j,k)+ dym*dx*pp(k)
124          jo_grad_x % p(i,j+1,k)   = jo_grad_x % p(i,j+1,k)+ dy*dxm*pp(k)
125          jo_grad_x % p(i+1,j+1,k) = jo_grad_x % p(i+1,j+1,k)+dy*dx*pp(k)
126       end do
127    end do                  
129    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_bogus_adj")
131 end subroutine da_transform_xtoy_bogus_adj