da_transfer_xatokma.inc

References to this file elsewhere.
1 subroutine da_transfer_xatokma(grid)
2 
3    !---------------------------------------------------------------------------
4    !  Purpose: Convert analysis increments into KMA increments 
5    !---------------------------------------------------------------------------
6 
7    implicit none
8    
9    type(domain),    intent(inout) :: grid
10 
11    integer :: i, j, k
12    integer :: is, ie, js, je, ks, ke
13    real    :: PU, PD, coeff
14 
15    if (trace_use) call da_trace_entry("da_transfer_xatokma")
16 
17    !------------------------------------------------------------------------
18    ! Set array range indices for processor subdomain.
19    !------------------------------------------------------------------------
20 
21    is = grid%xp % its
22    ie = grid%xp % ite
23    js = grid%xp % jts
24    je = grid%xp % jte
25    ks = grid%xp % kts
26    ke = grid%xp % kte
27 
28    !---------------------------------------------------------------------------
29    ! Add increment to the original guess and update xb and "grid"
30    !---------------------------------------------------------------------------
31 
32    do j=js,je
33       do i=is,ie
34          grid%xb%w(i,j,ke+1)=  grid%xb%w(i,j,ke+1) + grid%xa%w(i,j,ke+1)
35       end do
36       do i=is,ie
37          do k = ks, ke
38             grid%xb%u(i,j,k)   = grid%xa%u(i,j,k) + grid%xb%u(i,j,k)
39             grid%xb%v(i,j,k)   = grid%xa%v(i,j,k) + grid%xb%v(i,j,k)
40             grid%xb%t(i,j,k)   = grid%xa%t(i,j,k) + grid%xb%t(i,j,k)
41             grid%xb%w(i,j,k)   = grid%xa%w(i,j,k) + grid%xb%w(i,j,k)
42             grid%xb%q(i,j,k)   = grid%xa%q(i,j,k) + grid%xb%q(i,j,k)
43             ! compute pressure increments at KMA full levels
44             ! Note: Psfc is in hPa in  P = A + B * Psfc 
45             if (k == ke) then
46                coeff = grid%xb%KMA_B(K)/ &
47                   (grid%xb%KMA_A(K)+grid%xb%KMA_B(K)*grid%xb%psfc(I,J)/100.0)
48             else
49                PU = grid%xb%KMA_A(K+1) + &
50                   grid%xb%KMA_B(K+1)*grid%xb%psfc(I,J)/100.0
51                PD = grid%xb%KMA_A(K ) + &
52                   grid%xb%KMA_B(K )*grid%xb%psfc(I,J)/100.0
53                coeff=grid%xb%KMA_B(K) * &
54                   1.0/(PD-PU)**2*(-PU*(LOG(PD)-LOG(PU))+PD-PU) &
55                   + grid%xb%KMA_B(K+1)* &
56                   1.0/(PD-PU)**2*(PD*(LOG(PD)-LOG(PU))-PD+PU)
57             end if
58             grid%xa%p(i,j,k) = grid%xa%psfc(i,j) * coeff
59             grid%xa%p(i,j,k) = grid%xb%psfc(i,j)*grid%xa%psfc(i,j)
60             grid%xb%p(i,j,k) = grid%xb%p(i,j,k) + grid%xa%p(I,J,k)
61          end do
62          grid%xb%psfc(i,j) = grid%xb%psfc(i,j) + grid%xa%psfc(i,j)
63       end do
64    end do
65 
66    if (write_increments) call da_write_kma_increments(grid%xp, grid%xb, grid%xa)
67 
68    do j=js,je
69       do i=is,ie
70         grid%em_w_2(i,j,ke+1)=  grid%em_w_2(i,j,ke+1) + grid%xa%w(i,j,ke+1)
71         grid%psfc(i,j) = grid%psfc(i,j) + grid%xa%psfc(i,j)
72       end do
73    end do
74 
75    do k=ks,ke
76       do j=js,je
77          do i=is,ie
78             grid%em_u_2(i,j,k) = grid%em_u_2(i,j,k) + grid%xa%u(i,j,k)
79             grid%em_v_2(i,j,k) = grid%em_v_2(i,j,k) + grid%xa%v(i,j,k)
80             grid%em_w_2(i,j,k) = grid%em_w_2(i,j,k) + grid%xa%w(i,j,k)
81             grid%moist(i,j,k,P_QV) = grid%moist(i,j,k,P_QV) + grid%xa%q(i,j,k)
82          end do
83       end do
84    end do
85 
86    if (trace_use) call da_trace_exit("da_transfer_xatokma")
87 
88 end subroutine da_transfer_xatokma
89 
90