da_uv_to_divergence_adj.inc

References to this file elsewhere.
1 subroutine da_uv_to_divergence_adj(xp, xb, u, v, div)
2 
3    !---------------------------------------------------------------------------
4    !  Purpose: Adjoint of the subroutine da_uv_to_divergence
5    !                        d   U      d   V
6    !           Div = m^2 *[---(---) + ---(---) ] 
7    !                        dx  m      dy  M
8    !---------------------------------------------------------------------------
9 
10    implicit none
11 
12    type (xpose_type), intent(inout) :: xp      ! Dimensions and xpose buffers.
13 
14    type (xb_type), intent(in)           :: xb         ! First guess structure.
15    real, intent(out)  :: u(ims:ime,jms:jme,kms:kme)   ! u wind comp.
16    real, intent(out)  :: v(ims:ime,jms:jme,kms:kme)   ! v wind comp.
17    real, intent(inout):: div(ims:ime,jms:jme,kms:kme) ! Divergence.
18 
19    integer            :: i, j, k                      ! Loop counters.
20    integer            :: is, ie                       ! 1st dim. end points.
21    integer            :: js, je                       ! 2nd dim. end points.
22    real               :: one_third                    ! 1/3.
23 
24    real               :: coeff, inv_2ds
25    real               :: um(ims:ime,jms:jme)          ! Temp. storage of u/m.
26    real               :: vm(ims:ime,jms:jme)          ! Temp. storage of v/m
27 
28    !---------------------------------------------------------------------------
29    ! [1.0] Initialise:
30    !---------------------------------------------------------------------------
31 
32    if (.not. global) then
33       um(ims:ime,jms:jme) = 0.
34       vm(ims:ime,jms:jme) = 0.
35       inv_2ds = 0.5 / xb%ds
36    end if
37 
38    one_third = 1.0 / 3.0
39 
40    if (Testing_WRFVAR) then
41       is = its; ie = ite; js = jts; je = jte
42    else
43       is = its - 1; ie = ite + 1; js = jts - 1; je = jte + 1
44    end if
45 
46    if (.not. global .and. its == ids) is = ids+1
47    if (.not. global .and. ite == ide) ie = ide-1
48    if (jts == jds) js = jds+1; if (jte == jde) je = jde-1
49 
50    !---------------------------------------------------------------------------
51    ! [2.0] Calculate divergence:
52    !---------------------------------------------------------------------------
53 
54    if (global) then
55       call da_set_boundary_3d(div)
56    end if
57 
58    do k =kds, kde
59       if (.not. global) then
60          ! [2.2] Impose zero divergence gradient condition at boundaries:
61 
62          ! [2.2.4] Right boundaries:
63 
64          if (jte == jde) then
65             j = jte
66             do i = its, ite    ! This is different to original
67                div(i,j-1,k)=div(i,j-1,k)+div(i,j,k)*one_third*4.0
68                div(i,j-2,k)=div(i,j-2,k)-div(i,j,k)*one_third
69                div(i,j,k)=0.
70             end do
71          end if
72 
73          ! [2.2.3] Left boundaries:
74 
75          if (jts == jds) then
76             j = jts
77             do i = its, ite    ! This is different to original
78                div(i,j+1,k)=div(i,j+1,k)+div(i,j,k)*one_third*4.0
79                div(i,j+2,k)=div(i,j+2,k)-div(i,j,k)*one_third
80                div(i,j,k)=0.
81             end do
82          end if
83 
84          ! [2.2.2] Top boundaries:
85 
86          if (ite == ide) then
87             i = ite
88             do j = jts, jte
89                div(i-1,j,k)=div(i-1,j,k)+div(i,j,k)*one_third*4.0
90                div(i-2,j,k)=div(i-2,j,k)-div(i,j,k)*one_third
91                div(i,j,k)=0.
92             end do
93          end if
94 
95          ! [2.2.1] Bottom boundaries:
96 
97          if (its == ids) then
98             i = its 
99             do j = jts, jte
100                div(i+1,j,k)=div(i+1,j,k)+div(i,j,k)*one_third*4.0
101                div(i+2,j,k)=div(i+2,j,k)-div(i,j,k)*one_third
102                div(i,j,k)=0.
103             end do
104          end if
105 
106          ! [2.1] Compute fd divergence at interior points:
107          ! Computation to check for edge of domain:
108          ! This is only for adjoint, as we have to cross the processor boundary
109          ! to get the contribution.
110 
111          if (.not. Testing_WRFVAR) then
112             xp%vxy(its:ite, jts:jte) = div(its:ite, jts:jte, k)
113             call wrf_dm_halo(xp%domdesc,xp%comms,xp%halo_id12)
114 
115             div(is, js:je, k) = xp%vxy(is, js:je)
116             div(ie, js:je, k) = xp%vxy(ie, js:je)
117             div(is:ie, js, k) = xp%vxy(is:ie, js)
118             div(is:ie, je, k) = xp%vxy(is:ie, je)
119          end if
120 
121          do j = js, je
122             do i = is, ie
123               coeff = xb%map_factor(i,j) * xb%map_factor(i,j) * inv_2ds
124               um(i+1,j)=um(i+1,j)+div(i,j,k)*coeff
125               um(i-1,j)=um(i-1,j)-div(i,j,k)*coeff
126               vm(i,j+1)=vm(i,j+1)+div(i,j,k)*coeff
127               vm(i,j-1)=vm(i,j-1)-div(i,j,k)*coeff
128            end do
129          end do
130 
131          u(is-1:ie+1,js-1:je+1,k) = u(is-1:ie+1,js-1:je+1,k) &
132             +um(is-1:ie+1,js-1:je+1) / xb%map_factor(is-1:ie+1,js-1:je+1)
133          v(is-1:ie+1,js-1:je+1,k) = v(is-1:ie+1,js-1:je+1,k) &
134             +vm(is-1:ie+1,js-1:je+1) / xb%map_factor(is-1:ie+1,js-1:je+1)
135 
136          um(:,:)=0.
137          vm(:,:)=0.
138 
139       else
140 
141          !-------------------------------------------------------------------------
142          ! [2.1] Compute fd divergence at interior points:
143          !-------------------------------------------------------------------------
144 
145          do j = je, js, -1
146             do i = ie, is, -1  
147                u(i+1,j,k) = u(i+1,j,k) + xb%coefx(i,j) * div(i,j,k) 
148                u(i-1,j,k) = u(i-1,j,k) - xb%coefx(i,j) * div(i,j,k) 
149                v(i,j+1,k) = v(i,j+1,k) + xb%coefy(i,j) * div(i,j,k) 
150                v(i,j-1,k) = v(i,j-1,k) - xb%coefy(i,j) * div(i,j,k) 
151             end do
152          end do
153       end if
154    end do
155 
156    div = 0.0
157 
158 end subroutine da_uv_to_divergence_adj
159 
160