updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_vtox_transforms / da_calc_flow_dependence_xa.inc
blobc9ea003958192bc9ed33563771b53169f2764942
1 subroutine da_calc_flow_dependence_xa (xa_ens, ne,  ep, vp, nobwin)
3    !-----------------------------------------------------------------------
4    ! Purpose: calculate flow-dependent increments in model space (xa)
5    !          for a certain sub-window
6    !-----------------------------------------------------------------------
8    implicit none
10    type (x_subtype), intent(inout)      :: xa_ens
11    integer, intent(in)                  :: ne  ! Ensemble size.
12    type (ep_type), intent(in)           :: ep  ! Ensemble perturbations.
13    type (vp_type), intent(in)           :: vp  ! CV on grid structure.
14    integer, intent(in), optional        :: nobwin
16    integer                              :: iobwin, ii
17    integer                              :: i, j, k, n  ! Loop counters.
18    real                                 :: alpha       ! Local alpha copy.
20    if (trace_use) call da_trace_entry("da_calc_flow_dependence_xa")
22    iobwin = 1
23    if ( present(nobwin) ) then
24       iobwin = nobwin
25    end if
27    do n = 1, ne
28       ii = (iobwin-1)*ensdim_alpha + n
30       do k = kts, kte
31          do j = jts, jte
32             do i = its, ite
34                alpha = vp % alpha(i,j,k,n)
36                ! u:
37                xa_ens % u(i,j,k) = xa_ens % u(i,j,k) + alpha * ep % v1(i,j,k,ii) ! v1 = u
39                ! v:
40                xa_ens % v(i,j,k) = xa_ens % v(i,j,k) + alpha * ep % v2(i,j,k,ii) ! v2 = v
42                ! t:
43                xa_ens % t(i,j,k) = xa_ens % t(i,j,k) + alpha * ep % v3(i,j,k,ii) ! v3 = t
45                ! q:
46                xa_ens % q(i,j,k) = xa_ens % q(i,j,k) + alpha * ep % v4(i,j,k,ii) ! v4 = q
48             end do !i loop
49          end do !j loop
50       end do !k loop
52       xa_ens % psfc(its:ite,jts:jte) = xa_ens % psfc(its:ite,jts:jte) +  &
53                                        vp % alpha(its:ite,jts:jte,1,n) * &
54                                        ep % v5(its:ite,jts:jte,1,ii)       ! v5 = ps
56       if ( alpha_hydrometeors ) then
57          do k = kts, kte
58             do j = jts, jte
59                do i = its, ite
60                   alpha = vp % alpha(i,j,k,n)
61                   xa_ens%qcw(i,j,k) = xa_ens % qcw(i,j,k) + alpha * ep%cw(i,j,k,ii)
62                   xa_ens%qrn(i,j,k) = xa_ens % qrn(i,j,k) + alpha * ep%rn(i,j,k,ii)
63                   xa_ens%qci(i,j,k) = xa_ens % qci(i,j,k) + alpha * ep%ci(i,j,k,ii)
64                   xa_ens%qsn(i,j,k) = xa_ens % qsn(i,j,k) + alpha * ep%sn(i,j,k,ii)
65                   xa_ens%qgr(i,j,k) = xa_ens % qgr(i,j,k) + alpha * ep%gr(i,j,k,ii)
66                end do !i loop
67             end do !j loop
68          end do !k loop
69       end if ! alpha_hydrometeors
71    end do !n loop
73    if (trace_use) call da_trace_exit("da_calc_flow_dependence_xa")
75 end subroutine da_calc_flow_dependence_xa