xref: /phasta/phSolver/common/setSuction_Duct3.f90 (revision 712d3df0b59ebebaaeaea358162c8d2c043c6e08)
110167291SKenneth E. Jansenc========================================================================
210167291SKenneth E. Jansenc Dynamically controll the suction level based on main contraction mdot,
310167291SKenneth E. Jansenc called by itrdrv.f
410167291SKenneth E. Jansenc========================================================================
510167291SKenneth E. Jansen        subroutine setSuction_Duct3(x,BC,y, ilwork)
610167291SKenneth E. Jansen
710167291SKenneth E. Jansen        use wallData ! wnorm
810167291SKenneth E. Jansen        include "common.h"
910167291SKenneth E. Jansen        include "mpif.h"
1010167291SKenneth E. Jansen        include "auxmpi.h"
1110167291SKenneth E. Jansen
1210167291SKenneth E. Jansen        integer i, j, k, nn, nL
1310167291SKenneth E. Jansen        integer ilwork(nlwork)
1410167291SKenneth E. Jansen        integer nSuctionNode
1510167291SKenneth E. Jansen        integer, allocatable :: suctionNodeMap(:)
1610167291SKenneth E. Jansen!        logical inPatch
1710167291SKenneth E. Jansen        real*8 xcoor,ycoor,zcoor
1810167291SKenneth E. Jansen        real*8 BC(nshg,ndofBC)
1910167291SKenneth E. Jansen        real*8 BC3(numnp, 5)                      !used for hack to get suction right on part boundaries
2010167291SKenneth E. Jansen        real*8 x(nshg,nsd), y(nshg,ndof)
2110167291SKenneth E. Jansen!       real*8 xmin_t, xmax_t, xmin_b, xmax_b     !ramp geometry controls
2210167291SKenneth E. Jansen!       real*8 xmin_su, xmax_su, xmin_sl, xmax_sl !side wall patch limits
2310167291SKenneth E. Jansen!       real*8 t_upper                            !upper thickness of suction area
2410167291SKenneth E. Jansen!       real*8 ymin_su, ymax_sl
2510167291SKenneth E. Jansen        real*8 nvel_s, nvel_t, nvel_b             !normal velocity on the side walls, top, and bottom
2610167291SKenneth E. Jansen
2710167291SKenneth E. Jansen!        real*8 delX, delY, X1, Y1     !Variables used to map upper surface or trumpet
2810167291SKenneth E. Jansen!        real*8 Xline(30,2)
2910167291SKenneth E. Jansen!        integer suctionIDs(4)
3010167291SKenneth E. Jansen
3110167291SKenneth E. Jansen!        real*8 wallBCout(nshg, 6)
3210167291SKenneth E. Jansen
3310167291SKenneth E. Jansen        if(myrank.eq.0)write(*,*)'Setting side suctions'
3410167291SKenneth E. Jansen
3510167291SKenneth E. Jansen        nvel_b     = suctionVbottom      !normal velocity at the bottom patch
3610167291SKenneth E. Jansen        nvel_sl    = suctionVside_lower  !normal velocity of lower side wall patches
3710167291SKenneth E. Jansen        nvel_su    = suctionVside_upper  !normal velocity of upper side wall patches
3810167291SKenneth E. Jansen        nvel_t     = suctionVtop         !normal velocity at the top patch
3910167291SKenneth E. Jansen
4010167291SKenneth E. Jansen        allocate(suctionNodeMap(nshg))
4110167291SKenneth E. Jansen        call sfID2np(isetSuctionID_Duct, nSuctionNode, suctionNodeMap)  !to figure out whether it's the top, bottom, or side surface.
4210167291SKenneth E. Jansen
4310167291SKenneth E. Jansen        smallNum=1.0e-5
4410167291SKenneth E. Jansen        xsideBotMax=0.0428625+smallNum
4510167291SKenneth E. Jansen        xsideTopMax=0.130175+smallNum
4610167291SKenneth E. Jansen
4710167291SKenneth E. Jansen        do i = 1, nSuctionNode
4810167291SKenneth E. Jansen          nn = suctionNodeMap(i) ! map face node index to global node
4910167291SKenneth E. Jansen          xcoor = x(nn,1)
5010167291SKenneth E. Jansen          ycoor = x(nn,2)
5110167291SKenneth E. Jansen          !zcoor = x(nn,3)
5210167291SKenneth E. Jansen
5310167291SKenneth E. Jansen          !wnormx = wnorm(nn, 1)
5410167291SKenneth E. Jansen          wnormy = wnorm(nn, 2)
5510167291SKenneth E. Jansen          wnormz = wnorm(nn, 3)
5610167291SKenneth E. Jansen
5710167291SKenneth E. Jansen          !Test if the point lies in one of the suction patches
5810167291SKenneth E. Jansen
5910167291SKenneth E. Jansen          if((abs(nvel_sl) >  0  .or. abs(nvel_su) > 0) .and.  !if either upper or lower side walls are on...
6010167291SKenneth E. Jansen     &       (abs(wnormz) >=  0.999)) then                     ! and the node normal is in the z direction
6110167291SKenneth E. Jansen            if((ycoor > 0) .and. (xcoor .le. xsideTopMax)) then !hack to prevent edges from turning on.
6210167291SKenneth E. Jansen              BC(nn, 3) = nvel_su*wnorm(nn, 1) ! set xVel
6310167291SKenneth E. Jansen              BC(nn, 4) = nvel_su*wnorm(nn, 2) ! set yVel
6410167291SKenneth E. Jansen              BC(nn, 5) = nvel_su*wnorm(nn, 3) ! set zVel
6510167291SKenneth E. Jansen            elseif (xcoor .le. xsideBotMax) then
6610167291SKenneth E. Jansen              BC(nn, 3) = nvel_sl*wnorm(nn, 1) ! set xVel
6710167291SKenneth E. Jansen              BC(nn, 4) = nvel_sl*wnorm(nn, 2) ! set yVel
6810167291SKenneth E. Jansen              BC(nn, 5) = nvel_sl*wnorm(nn, 3) ! set zVel
6910167291SKenneth E. Jansen            endif
7010167291SKenneth E. Jansen          else if(abs(wnormz) > 0.1) then      ! side walls should have a larger znormal than this - it's ambiguous what to do so do nothing
7110167291SKenneth E. Jansen            cycle
7210167291SKenneth E. Jansen          else if(abs(nvel_t) > 0      .and. !top suction patch
7310167291SKenneth E. Jansen     &                 wnormy > 0.001) then  !face is pointing up (in the +y direction)
7410167291SKenneth E. Jansen            BC(nn, 3) = nvel_t*wnorm(nn, 1) ! set xVel using wall normal
7510167291SKenneth E. Jansen            BC(nn, 4) = nvel_t*wnorm(nn, 2) ! set yVel
7610167291SKenneth E. Jansen            BC(nn, 5) = nvel_t*wnorm(nn, 3) ! set zVel
7710167291SKenneth E. Jansen
7810167291SKenneth E. Jansen          else if(abs(nvel_b) > 0       .and. !bottom suction patches
7910167291SKenneth E. Jansen     &                wnormy  < -0.001) then  !face is pointing down (in the -y direction)
8010167291SKenneth E. Jansen            BC(nn, 3) = nvel_b*wnorm(nn, 1) ! set xVel using wall normal
8110167291SKenneth E. Jansen            BC(nn, 4) = nvel_b*wnorm(nn, 2) ! set yVel
8210167291SKenneth E. Jansen            BC(nn, 5) = nvel_b*wnorm(nn, 3) ! set zVel
8310167291SKenneth E. Jansen          endif
8410167291SKenneth E. Jansen
8510167291SKenneth E. Jansen        enddo
8610167291SKenneth E. Jansen
87*8636783dSKenneth E. Jansen        BC3=zero
8810167291SKenneth E. Jansen        if(numpe.gt.1) then
8910167291SKenneth E. Jansen            BC3(:,1:3)=BC(:,3:5)
9010167291SKenneth E. Jansen            do i=1,nshg
9110167291SKenneth E. Jansen                bc3mag=BC3(i,1)**2+BC3(i,2)**2+BC3(i,3)**2
9210167291SKenneth E. Jansen                if(bc3mag.gt.0) BC3(i,4)=one
9310167291SKenneth E. Jansen            enddo
9410167291SKenneth E. Jansen            call commu(BC3,ilwork,5,'in ')
9510167291SKenneth E. Jansen! This accumulated and failed                 call commu(BC(:,3:5),ilwork,3,'in ')
9610167291SKenneth E. Jansen            do i=1,nshg
9710167291SKenneth E. Jansen                bcmag=BC(i,3)**2+BC(i,4)**2+BC(i,5)**2
9810167291SKenneth E. Jansen                if((bcmag.eq.0).and.(BC3(i,4).ne.0)) then !node is slave and has not been correctly set.
9910167291SKenneth E. Jansen                     BC(i,3:5)=BC3(i,1:3)/BC3(i,4)  !division by BC3 is necessary to account for more than two blocks sharing the same node
10010167291SKenneth E. Jansen                endif
10110167291SKenneth E. Jansen            enddo
10210167291SKenneth E. Jansen       endif
10310167291SKenneth E. Jansen
10410167291SKenneth E. Jansen        !Debugging loop
10510167291SKenneth E. Jansen!        if(nSuctionNode.gt.0) then
10610167291SKenneth E. Jansen!          icount=0
10710167291SKenneth E. Jansen!          do i=1,nSuctionNode
10810167291SKenneth E. Jansen!             nn=suctionNodeMap(i)
10910167291SKenneth E. Jansen!             dz=dabs(dabs(x(nn,3))-0.0762)
11010167291SKenneth E. Jansen!             if(dz.lt.0.0001) then
11110167291SKenneth E. Jansen!               icount=icount+1
11210167291SKenneth E. Jansen!               if(icount.eq.1)
11310167291SKenneth E. Jansen!     & write(myrank+1000,*) 'writing wnorm and BC setSuctionNG3'
11410167291SKenneth E. Jansen!               write(myrank+1000,1234) nn,i,
11510167291SKenneth E. Jansen!     &           x(nn,1), x(nn,2), x(nn,3),
11610167291SKenneth E. Jansen!     &           wnorm(nn,1),wnorm(nn,2),wnorm(nn,3),
11710167291SKenneth E. Jansen!     &           BC(nn,3),BC(nn,4),BC(nn,5)
11810167291SKenneth E. Jansen!             endif
11910167291SKenneth E. Jansen!          enddo
12010167291SKenneth E. Jansen!         if(icount.gt.0) close(myrank+1000)
12110167291SKenneth E. Jansen!        endif
12210167291SKenneth E. Jansen!1234           format(i5,2x,i5,9(2x,e14.7))
12310167291SKenneth E. Jansen
12410167291SKenneth E. Jansen      end
12510167291SKenneth E. Jansen
12610167291SKenneth E. Jansen
12710167291SKenneth E. Jansen
12810167291SKenneth E. Jansen
12910167291SKenneth E. Jansen
13010167291SKenneth E. Jansen
131