da_transform_xtoy_pseudo_adj.inc
References to this file elsewhere.
1 subroutine da_transform_xtoy_pseudo_adj(iv, jo_grad_y, jo_grad_x)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 !-----------------------------------------------------------------------
6
7 implicit none
8
9 type (iv_type), intent(in) :: iv ! obs. inc vector (o-b).
10 type (y_type) , intent(in) :: jo_grad_y ! grad_y(jo)
11 type (x_type) , intent(inout) :: jo_grad_x ! grad_x(jo)
12
13 integer :: n ! Loop counter.
14
15 real, allocatable :: u(:,:)
16 real, allocatable :: v(:,:)
17 real, allocatable :: q(:,:)
18 real, allocatable :: p(:,:)
19 real, allocatable :: t(:,:)
20
21 if (trace_use_dull) call da_trace_entry("da_transform_xtoy_pseudo_adj")
22
23 allocate (u(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
24 allocate (v(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
25 allocate (q(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
26 allocate (p(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
27 allocate (t(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
28
29 do n=iv%info(pseudo)%n1,iv%info(pseudo)%n2
30 u(1,n) = jo_grad_y%pseudo(n)%u
31 v(1,n) = jo_grad_y%pseudo(n)%v
32 q(1,n) = jo_grad_y%pseudo(n)%q
33 p(1,n) = jo_grad_y%pseudo(n)%p
34 t(1,n) = jo_grad_y%pseudo(n)%t
35 end do
36
37 call da_interp_lin_3d_adj(jo_grad_x%u, iv%info(pseudo), u)
38 call da_interp_lin_3d_adj(jo_grad_x%v, iv%info(pseudo), v)
39 call da_interp_lin_3d_adj(jo_grad_x%q, iv%info(pseudo), q)
40 call da_interp_lin_3d_adj(jo_grad_x%p, iv%info(pseudo), p)
41 call da_interp_lin_3d_adj(jo_grad_x%t, iv%info(pseudo), t)
42
43 deallocate (u)
44 deallocate (v)
45 deallocate (q)
46 deallocate (p)
47 deallocate (t)
48
49 if (trace_use_dull) call da_trace_exit("da_transform_xtoy_pseudo_adj")
50
51 end subroutine da_transform_xtoy_pseudo_adj
52
53