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
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 if (jts == jds) then
35 do k = kts,kte
36 tmp_u(k) = tmpvar*sum(-u(its:ite,jts,k)*cos_xls(its:ite) &
37 +v(its:ite,jts,k)*sin_xls(its:ite))
38 tmp_v(k) = tmpvar*sum(-u(its:ite,jts,k)*sin_xls(its:ite) &
39 -v(its:ite,jts,k)*cos_xls(its:ite))
40 end do
41 end if
42
43 call wrf_dm_sum_reals(tmp_u(:), tmpu(:))
44 call wrf_dm_sum_reals(tmp_v(:), tmpv(:))
45
46 if (jts == jds) then
47 do k = kts,kte
48 u(its:ite,jts+1,k) = u(its:ite,jts+1,k) -tmpu(k)*cos_xls(its:ite) &
49 - tmpv(k)*sin_xls(its:ite)
50 v(its:ite,jts+1,k) = v(its:ite,jts+1,k) +tmpu(k)*sin_xls(its:ite) &
51 - tmpv(k)*cos_xls(its:ite)
52 u(its:ite,jts,k) = 0.0
53 v(its:ite,jts,k) = 0.0
54 end do
55 end if
56
57 tmp_u(:) =0.0
58 tmp_v(:) =0.0
59
60 if (jte == jde) then
61 do k = kts,kte
62 tmp_u(k) = tmpvar*sum(-u(its:ite,jte,k)*cos_xle(its:ite) &
63 -v(its:ite,jte,k)*sin_xle(its:ite))
64 tmp_v(k) = tmpvar*sum( u(its:ite,jte,k)*sin_xle(its:ite) &
65 -v(its:ite,jte,k)*cos_xle(its:ite))
66 end do
67 end if
68
69 call wrf_dm_sum_reals(tmp_u(:), tmpu(:))
70 call wrf_dm_sum_reals(tmp_v(:), tmpv(:))
71
72 if (jte == jde) then
73 do k = kts,kte
74 u(its:ite,jte-1,k) = u(its:ite,jte-1,k) -tmpu(k)*cos_xle(its:ite) &
75 + tmpv(k)*sin_xle(its:ite)
76 v(its:ite,jte-1,k) = v(its:ite,jte-1,k) -tmpu(k)*sin_xle(its:ite) &
77 - tmpv(k)*cos_xle(its:ite)
78 u(its:ite,jte,k) = 0.0
79 v(its:ite,jte,k) = 0.0
80 end do
81 end if
82
83 if (trace_use) call da_trace_exit("da_get_avpoles")
84
85 end subroutine da_get_avpoles
86
87