da_get_avpoles.inc

References to this file elsewhere.
1 subroutine da_get_avpoles(u,v,       &
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 Polar winds 
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) :: u(ims:ime,jms:jme,kms:kme)   ! u wind comp.
16    real,    intent(inout) :: v(ims:ime,jms:jme,kms:kme)   ! v wind comp.
17  
18    real                   :: tmpvar                                         
19    real                   :: tmpu(kts:kte)
20    real                   :: tmp_u(kts:kte)
21    real                   :: tmpv(kts:kte)
22    real                   :: tmp_v(kts:kte)
23    integer                :: k, n
24 
25    if (trace_use) call da_trace_entry("da_get_avpoles")
26 
27    tmpvar      = 1.0/real(ide-ids+1)
28 
29    ! cos_xls etc in da_control, calculated in da_setup_firstguess
30 
31    tmp_u(:) =0.0
32    tmp_v(:) =0.0
33 
34    n = kte-kts+1
35 
36    if (jts == jds) then 
37       do k = kts,kte
38          tmp_u(k) = tmpvar*sum(-u(its:ite,jts,k)*cos_xls(its:ite) & 
39                             +v(its:ite,jts,k)*sin_xls(its:ite))
40          tmp_v(k) = tmpvar*sum(-u(its:ite,jts,k)*sin_xls(its:ite) & 
41                             -v(its:ite,jts,k)*cos_xls(its:ite))
42       end do
43    end if
44 
45    call wrf_dm_sum_reals(tmp_u(:), tmpu(:))
46    call wrf_dm_sum_reals(tmp_v(:), tmpv(:))
47 
48    if (jts == jds) then 
49       do k = kts,kte
50          u(its:ite,jts+1,k) = u(its:ite,jts+1,k) -tmpu(k)*cos_xls(its:ite) &
51                               - tmpv(k)*sin_xls(its:ite)
52          v(its:ite,jts+1,k) = v(its:ite,jts+1,k) +tmpu(k)*sin_xls(its:ite) &
53                               - tmpv(k)*cos_xls(its:ite)
54          u(its:ite,jts,k) = 0.0
55          v(its:ite,jts,k) = 0.0
56       end do
57    end if
58 
59    tmp_u(:) =0.0
60    tmp_v(:) =0.0
61 
62    if (jte == jde) then 
63       do k = kts,kte
64          tmp_u(k) = tmpvar*sum(-u(its:ite,jte,k)*cos_xle(its:ite) & 
65                             -v(its:ite,jte,k)*sin_xle(its:ite))
66          tmp_v(k) = tmpvar*sum( u(its:ite,jte,k)*sin_xle(its:ite) & 
67                             -v(its:ite,jte,k)*cos_xle(its:ite))
68       end do
69    end if
70 
71    call wrf_dm_sum_reals(tmp_u(:), tmpu(:))
72    call wrf_dm_sum_reals(tmp_v(:), tmpv(:))
73 
74    if (jte == jde) then 
75       do k = kts,kte
76          u(its:ite,jte-1,k) = u(its:ite,jte-1,k) -tmpu(k)*cos_xle(its:ite) &
77                               + tmpv(k)*sin_xle(its:ite)
78          v(its:ite,jte-1,k) = v(its:ite,jte-1,k) -tmpu(k)*sin_xle(its:ite) &
79                               - tmpv(k)*cos_xle(its:ite)
80          u(its:ite,jte,k) = 0.0
81          v(its:ite,jte,k) = 0.0
82       end do
83    end if
84 
85    if (trace_use) call da_trace_exit("da_get_avpoles")
86 
87 end subroutine da_get_avpoles
88 
89