da_get_aspoles.inc
References to this file elsewhere.
1 subroutine da_get_aspoles(x, &
2 ids, ide, jds, jde, &
3 ims, ime, jms, jme, kms, kme, &
4 its, ite, jts, jte, kts, kte )
5
6 !---------------------------------------------------------------------------
7 ! Purpose: Treatment for Adjoint of Scalar field at Poles
8 !---------------------------------------------------------------------------
9
10 implicit none
11
12 integer, intent(in) :: ids, ide, jds, jde
13 integer, intent(in) :: ims, ime, jms, jme, kms, kme
14 integer, intent(in) :: its, ite, jts, jte, kts, kte
15 real, intent(inout) :: x(ims:ime,jms:jme,kms:kme)
16
17 integer :: k
18 real :: tmpvar
19 real :: tmps(kts:kte)
20 real :: tmp_s(kts:kte)
21
22 tmpvar = 1.0/real(ide-ids+1)
23
24 if (trace_use) call da_trace_entry("da_get_aspoles")
25
26 tmp_s(:) = 0.0
27
28 if (jts == jds) then
29 do k = kts, kte
30 tmp_s(k) = tmpvar*sum(x(its:ite,jts,k))
31 end do
32 end if
33
34 call wrf_dm_sum_reals(tmp_s(:),tmps(:))
35
36 if (jts == jds) then
37 do k = kts, kte
38 x(its:ite,jts+1,k) = x(its:ite,jts+1,k) + tmps(k)
39 x(its:ite,jts,k) = 0.0
40 end do
41 end if
42
43 tmp_s(:) = 0.0
44
45 if (jte == jde) then
46 do k = kts, kte
47 tmp_s(k) = tmpvar*sum(x(its:ite,jte,k))
48 end do
49 end if
50
51 call wrf_dm_sum_reals(tmp_s(:),tmps(:))
52
53 if (jte == jde) then
54 do k = kts, kte
55 x(its:ite,jte-1,k) = x(its:ite,jte-1,k) + tmps(k)
56 x(its:ite,jte,k) = 0.0
57 end do
58 end if
59
60 if (trace_use) call da_trace_exit("da_get_aspoles")
61
62 end subroutine da_get_aspoles
63
64