da_interp_lin_3d_adj.inc
References to this file elsewhere.
1 subroutine da_interp_lin_3d_adj(fm3d, info, fo3d)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 !-----------------------------------------------------------------------
6
7 implicit none
8
9 real, intent(inout) :: fm3d(ims:ime,jms:jme,kms:kme)
10 type(infa_type), intent(in) :: info
11 real, intent(in) :: fo3d(1:info%max_lev,info%n1:info%n2)
12
13 integer :: n,k
14 real :: fmz(kms:kme)
15
16 if (trace_use) call da_trace_entry("da_interp_lin_3d_adj")
17
18 do n=info%n1,info%n2
19 fmz = 0.0
20 do k = 1, info%levels(n)
21 if (info%k(k,n) > 0) then
22 fmz(info%k(k,n)) = info%dzm(k,n) * fo3d(k,n) + fmz(info%k(k,n))
23 fmz(info%k(k,n)+1) = info%dz (k,n) * fo3d(k,n) + fmz(info%k(k,n)+1)
24 end if
25 end do
26
27 do k=kts,kte
28 fm3d(info%i(k,n) ,info%j(k,n),k) = info%dym(k,n)*info%dxm(k,n)*fmz(k) + fm3d(info%i(k,n) ,info%j(k,n) ,k)
29 fm3d(info%i(k,n)+1,info%j(k,n),k) = info%dym(k,n)*info%dx (k,n)*fmz(k) + fm3d(info%i(k,n)+1,info%j(k,n) ,k)
30 fm3d(info%i(k,n) ,info%j(k,n)+1,k) = info%dy (k,n)*info%dxm(k,n)*fmz(k) + fm3d(info%i(k,n) ,info%j(k,n)+1,k)
31 fm3d(info%i(k,n)+1,info%j(k,n)+1,k) = info%dy (k,n)*info%dx (k,n)*fmz(k) + fm3d(info%i(k,n)+1,info%j(k,n)+1,k)
32 end do
33 end do
34
35 if (trace_use) call da_trace_exit("da_interp_lin_3d_adj")
36
37 end subroutine da_interp_lin_3d_adj
38
39