da_get_vpoles.inc
References to this file elsewhere.
1 subroutine da_get_vpoles(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 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,tmp_u,tmpv,tmp_v
20 integer :: k
21
22 if (trace_use) call da_trace_entry("da_get_vpoles")
23
24 ! cos_xls etc in da_control, calculated in da_setup_firstguess
25
26 tmpvar = 1.0/real(ide-ids+1)
27
28 do k = kts,kte
29 tmp_u =0.0
30 tmp_v =0.0
31 tmpu = 0.0
32 tmpv = 0.0
33
34 if (jts == jds) then
35 tmp_u = tmpvar*sum(-u(its:ite,jts+1,k)*cos_xls(its:ite)&
36 +v(its:ite,jts+1,k)*sin_xls(its:ite))
37 tmp_v = tmpvar*sum(-u(its:ite,jts+1,k)*sin_xls(its:ite)&
38 -v(its:ite,jts+1,k)*cos_xls(its:ite))
39 end if
40
41 tmpu = wrf_dm_sum_real( tmp_u)
42 tmpv = wrf_dm_sum_real( tmp_v)
43
44 if (jts == jds) then
45 u(its:ite,jts,k) = -tmpu*cos_xls(its:ite) -tmpv*sin_xls(its:ite)
46 v(its:ite,jts,k) = tmpu*sin_xls(its:ite) -tmpv*cos_xls(its:ite)
47 end if
48
49 tmp_u =0.0
50 tmp_v =0.0
51 tmpu = 0.0
52 tmpv = 0.0
53
54 if (jte == jde) then
55 tmp_u = tmpvar*sum(-u(its:ite,jte-1,k)*cos_xle(its:ite)&
56 -v(its:ite,jte-1,k)*sin_xle(its:ite))
57 tmp_v = tmpvar*sum( u(its:ite,jte-1,k)*sin_xle(its:ite)&
58 -v(its:ite,jte-1,k)*cos_xle(its:ite))
59 end if
60 tmpu = wrf_dm_sum_real( tmp_u)
61 tmpv = wrf_dm_sum_real( tmp_v)
62 if (jte == jde) then
63 u(its:ite,jte,k) = -tmpu*cos_xle(its:ite) +tmpv*sin_xle(its:ite)
64 v(its:ite,jte,k) = -tmpu*sin_xle(its:ite) -tmpv*cos_xle(its:ite)
65 end if
66 end do
67
68 if (trace_use) call da_trace_exit("da_get_vpoles")
69
70 end subroutine da_get_vpoles
71
72