updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_vtox_transforms / da_get_spoles.inc
blob8ef31a7bdf82cd6d7fde10c61913c306ea5a75a7
1 subroutine da_get_spoles(x,              &
2           ids, ide, jds, jde, &
3           ims, ime, jms, jme, kms, kme,  &
4           its, ite, jts, jte, kts, kte  )
6    !---------------------------------------------------------------------------
7    ! Purpose: Treatment for Scalar field at Poles
8    !---------------------------------------------------------------------------
10    implicit none
11    
13    integer, intent(in)    :: ids, ide, jds, jde
14    integer, intent(in)    :: ims, ime, jms, jme, kms, kme
15    integer, intent(in)    :: its, ite, jts, jte, kts, kte
16    real,    intent(inout) :: x(ims:ime,jms:jme,kms:kme)   
18    integer                :: k
19    real                   :: tmpvar,tmps,tmp_s
21    if (trace_use) call da_trace_entry("da_get_spoles")
23    tmpvar      = 1.0/real(ide-ids+1)
25    do k = kts, kte
26       tmps =0.0  ; tmp_s =0.0
27       if (jts == jds)  tmp_s = tmpvar*sum(x(its:ite,jts+1,k))
29       tmps = wrf_dm_sum_real( tmp_s)
30       if (jts == jds) x(its:ite,jts,k) = tmps
32       tmps =0.0  ; tmp_s =0.0
33       if (jte == jde) tmp_s = tmpvar*sum(x(its:ite,jte-1,k))
35       tmps = wrf_dm_sum_real( tmp_s)
36       if (jte == jde) x(its:ite,jte,k) = tmps
37    end do
39    if (trace_use) call da_trace_exit("da_get_spoles")
41 end subroutine da_get_spoles