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