da_interp_lin_3d_adj.inc
References to this file elsewhere.
1 subroutine da_interp_lin_3d_adj(fm3d, xp, &
2 i, j, dx, dy, dxm, dym, &
3 fo3d, ml, zk, nl)
4
5 !-----------------------------------------------------------------------
6 ! Purpose: TBD
7 !-----------------------------------------------------------------------
8
9 implicit none
10
11 type (xpose_type), intent(in) :: xp ! Dimensions and xpose buffers.
12 integer, intent(in) :: i, j
13 real, intent(in) :: dx, dxm, dy, dym
14 integer, intent(in) :: nl,ml
15 real, dimension(ml), intent(in) :: fo3d
16 real, dimension(nl), intent(in) :: zk
17 real, dimension(xp%ims:xp%ime,xp%jms:xp%jme,xp%kms:xp%kme), &
18 intent(inout) :: fm3d ! Input/Output variable
19
20 integer :: kk, k
21 real :: dz, dzm
22 real, dimension(xp%kms:xp%kme) :: fmz
23
24 if (trace_use_frequent) call da_trace_entry("da_interp_lin_3d_adj")
25
26 fmz = 0.0
27
28 do kk = 1, nl
29 if (zk(kk) > 0.0) then
30 call da_togrid(zk(kk), xp%kts, xp%kte, k, dz, dzm)
31 fmz(k) = dzm*fo3d(kk) + fmz(k)
32 fmz(k+1) = dz *fo3d(kk) + fmz(k+1)
33 end if
34 end do
35
36 fm3d(i ,j ,xp%kts:xp%kte) = dym*dxm*fmz(xp%kts:xp%kte) &
37 + fm3d(i ,j ,xp%kts:xp%kte)
38 fm3d(i+1,j ,xp%kts:xp%kte) = dym*dx *fmz(xp%kts:xp%kte) &
39 + fm3d(i+1,j ,xp%kts:xp%kte)
40 fm3d(i ,j+1,xp%kts:xp%kte) = dy *dxm*fmz(xp%kts:xp%kte) &
41 + fm3d(i ,j+1,xp%kts:xp%kte)
42 fm3d(i+1,j+1,xp%kts:xp%kte) = dy *dx *fmz(xp%kts:xp%kte) &
43 + fm3d(i+1,j+1,xp%kts:xp%kte)
44
45 if (trace_use_frequent) call da_trace_exit("da_interp_lin_3d_adj")
46
47 end subroutine da_interp_lin_3d_adj
48
49