da_transform_xtoy_bogus_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_bogus_adj(grid, iv, jo_grad_y, jo_grad_x)
2 
3    !---------------------------------------------------------------------
4    ! Purpose: the adjoint of bogus observation operators.
5    !---------------------------------------------------------------------
6 
7    implicit none
8 
9    type (domain),     intent(in)     :: grid
10    type (iv_type),    intent(in)     :: iv          ! obs. inc vector (o-b).
11    type (y_type) ,    intent(inout)  :: jo_grad_y   ! grad_y(jo)
12    type (x_type) ,    intent(inout)  :: jo_grad_x   ! grad_x(jo)
13 
14    integer :: n        ! Loop counter.
15    integer :: i, j, k  ! Index dimension.
16    integer :: num_levs
17    real    :: dx, dxm
18    real    :: dy, dym
19 
20    real    :: model_t(kms:kme)
21    real    :: tt(kms:kme)
22    real    :: model_q(kms:kme)
23    real    :: qq(kms:kme)
24    real    :: model_p_c(kms:kme)
25    real    :: pp(kms:kme)
26    real    :: psfc_loc,terr_loc,ppsfc   
27 
28    real, allocatable :: temp_u(:,:)
29    real, allocatable :: temp_v(:,:)
30    real, allocatable :: temp_t(:,:)
31    real, allocatable :: temp_q(:,:)           
32 
33    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_bogus_adj")
34 
35    allocate (temp_u(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
36    allocate (temp_v(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
37    allocate (temp_t(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
38    allocate (temp_q(iv%info(bogus)%max_lev,iv%info(bogus)%n1:iv%info(bogus)%n2))
39 
40    do n=iv%info(bogus)%n1,iv%info(bogus)%n2
41       temp_u(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%u
42       temp_v(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%v
43       temp_t(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%t
44       temp_q(1:size(jo_grad_y%bogus(n)%u),n)  = jo_grad_y%bogus(n)%q
45    end do
46 
47    ! [1.1] Adjoint feedback from Y to X for u and v:
48 
49    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(bogus), temp_u)
50    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(bogus), temp_v)
51    call da_interp_lin_3d_adj (jo_grad_x%t, iv%info(bogus), temp_t)
52    call da_interp_lin_3d_adj (jo_grad_x%q, iv%info(bogus), temp_q)
53    deallocate (temp_u)
54    deallocate (temp_v)
55    deallocate (temp_t)
56    deallocate (temp_q)
57 
58    do n= iv%info(bogus)%n1, iv%info(bogus)%n2
59       num_levs = iv%info(bogus)%levels(n)
60       if (num_levs < 1) cycle
61 
62       ! [1.1] Get cross pt. horizontal interpolation weights:
63 
64       i   = iv%info(bogus)%i(1,n)
65       dy  = iv%info(bogus)%dy(1,n)
66       dym = iv%info(bogus)%dym(1,n)
67       j   = iv%info(bogus)%j(1,n)
68       dx  = iv%info(bogus)%dx(1,n)
69       dxm = iv%info(bogus)%dxm(1,n)
70 
71       ! [1.2] Compute the feedback from SLP to t and q:
72 
73       ! 1.2.1 Background at the observation location:
74 
75       do k = kts, kte
76          model_t(k) = dym*(dxm*grid%xb%t(i,j,k) + dx*grid%xb%t(i+1,j,k)) &
77             + dy *(dxm*grid%xb%t(i,j+1,k) + dx*grid%xb%t(i+1,j+1,k))
78          model_q(k) = dym*(dxm*grid%xb%q(i,j,k) + dx*grid%xb%q(i+1,j,k)) &
79             + dy *(dxm*grid%xb%q(i,j+1,k) + dx*grid%xb%q(i+1,j+1,k))
80          model_p_c(k) = dym*(dxm*grid%xb%p(i,j,k) + dx*grid%xb%p(i+1,j,k)) &
81             + dy *(dxm*grid%xb%p(i,j+1,k) + dx*grid%xb%p(i+1,j+1,k))
82       end do
83 
84       terr_loc = dym*(dxm*grid%xb%terr(i,j)   + dx*grid%xb%terr(i+1,j)) &
85          + dy *(dxm*grid%xb%terr(i,j+1) + dx*grid%xb%terr(i+1,j+1))
86       psfc_loc = dym*(dxm*grid%xb%psfc(i,j)   + dx*grid%xb%psfc(i+1,j)) &
87          + dy *(dxm*grid%xb%psfc(i,j+1) + dx*grid%xb%psfc(i+1,j+1))
88 
89       ! 1.2.2 Compute the feedback from SLP to p, t, q, and psfc 
90       !       at the observed location:
91 
92       call da_tpq_to_slp_adj(model_t, model_q ,model_p_c, terr_loc, psfc_loc, &
93          tt, qq, pp, ppsfc, jo_grad_y%bogus(n)%slp)  
94 
95       ! 1.2.3 feedback from the observed location to grid space:
96 
97       ! 1.2.3.1 for Psfc
98 
99       jo_grad_x % psfc(i,j)     = jo_grad_x % psfc(i,j) + dym*dxm*ppsfc
100       jo_grad_x % psfc(i+1,j)   = jo_grad_x % psfc(i+1,j)+ dym*dx*ppsfc
101       jo_grad_x % psfc(i,j+1)   = jo_grad_x % psfc(i,j+1)+ dy*dxm*ppsfc
102       jo_grad_x % psfc(i+1,j+1) = jo_grad_x % psfc(i+1,j+1)+dy*dx*ppsfc
103 
104       ! 1.2.3.2 for t, q, p
105 
106       do k = kts, kte
107          jo_grad_x % t(i,j,k)     = jo_grad_x % t(i,j,k)+dym*dxm*tt(k)
108          jo_grad_x % t(i+1,j,k)   = jo_grad_x % t(i+1,j,k)+ dym*dx*tt(k)
109          jo_grad_x % t(i,j+1,k)   = jo_grad_x % t(i,j+1,k)+ dy*dxm*tt(k)
110          jo_grad_x % t(i+1,j+1,k) = jo_grad_x % t(i+1,j+1,k)+dy*dx*tt(k)
111          jo_grad_x % q(i,j,k)     = jo_grad_x % q(i,j,k)+dym*dxm*qq(k)
112          jo_grad_x % q(i+1,j,k)   = jo_grad_x % q(i+1,j,k)+ dym*dx*qq(k)
113          jo_grad_x % q(i,j+1,k)   = jo_grad_x % q(i,j+1,k)+ dy*dxm*qq(k)
114          jo_grad_x % q(i+1,j+1,k) = jo_grad_x % q(i+1,j+1,k)+dy*dx*qq(k)
115          jo_grad_x % p(i,j,k)     = jo_grad_x % p(i,j,k)+dym*dxm*pp(k)
116          jo_grad_x % p(i+1,j,k)   = jo_grad_x % p(i+1,j,k)+ dym*dx*pp(k)
117          jo_grad_x % p(i,j+1,k)   = jo_grad_x % p(i,j+1,k)+ dy*dxm*pp(k)
118          jo_grad_x % p(i+1,j+1,k) = jo_grad_x % p(i+1,j+1,k)+dy*dx*pp(k)
119       end do
120    end do                  
121 
122    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_bogus_adj")
123 
124 end subroutine da_transform_xtoy_bogus_adj
125 
126