xref: /phasta/phSolver/common/mod_wallData.f90 (revision 75f1c48cf62ef249f52322de5c51e7d33ddb5068)
110167291SKenneth E. Jansen      module wallData !was suctionDuct
210167291SKenneth E. Jansen        logical wnormCalced
310167291SKenneth E. Jansen        !integer nsuctionface
410167291SKenneth E. Jansen        !integer, allocatable :: suctionf(:)
510167291SKenneth E. Jansen
610167291SKenneth E. Jansen        !Threshold for determining whether a node is on a corner
710167291SKenneth E. Jansen        real*8 :: cornerCosThresh
810167291SKenneth E. Jansen
910167291SKenneth E. Jansen        !Wall normal array of size (nshg, 3)
1010167291SKenneth E. Jansen        real*8, allocatable :: wnorm(:,:)
11*75f1c48cSCameron Smith
1210167291SKenneth E. Jansen      contains
1310167291SKenneth E. Jansen
14*75f1c48cSCameron Smith        subroutine destroyWallData()
15*75f1c48cSCameron Smith          if ( allocated(wnorm) ) then
16*75f1c48cSCameron Smith            deallocate(wnorm)
17*75f1c48cSCameron Smith            wnormCalced = .false.
18*75f1c48cSCameron Smith          endif
19*75f1c48cSCameron Smith        end subroutine destroyWallData
20*75f1c48cSCameron Smith
2110167291SKenneth E. Jansen        subroutine findWallNorm(x,iBC,ilwork,iper)
2210167291SKenneth E. Jansen
2310167291SKenneth E. Jansen!         use wallData
2410167291SKenneth E. Jansen          use pointer_data
2510167291SKenneth E. Jansen          include "common.h"
2610167291SKenneth E. Jansen          include "mpif.h"
2710167291SKenneth E. Jansen
2810167291SKenneth E. Jansen          integer iBC(nshg)
2910167291SKenneth E. Jansen          real*8 x(numnp,nsd)
3010167291SKenneth E. Jansen          integer iper(nshg)
3110167291SKenneth E. Jansen          integer ilwork(nlwork)
3210167291SKenneth E. Jansen          integer, allocatable :: ienb(:)
3310167291SKenneth E. Jansen          real*8,  allocatable :: wNormTmp(:,:)
3410167291SKenneth E. Jansen          real*8 eNorm(nsd),en(nsd),wn(nsd),aedg(nsd),bedg(nsd)
3510167291SKenneth E. Jansen          real*8 mag
3610167291SKenneth E. Jansen!         integer isfID
3710167291SKenneth E. Jansen
3810167291SKenneth E. Jansen          real*8 wallBCout(nshg, 6)
3910167291SKenneth E. Jansen          real*8 BC(nshg,ndofBC)
4010167291SKenneth E. Jansen
4110167291SKenneth E. Jansen          cornerCosThresh = 0.75d00
4210167291SKenneth E. Jansen          !The number 0.75 is somewhat arbitrary. The thinking is that this
4310167291SKenneth E. Jansen          !code is also used for the normal on curved surfaces which will have
4410167291SKenneth E. Jansen          !a non unit dot product from elements surrounding a wall vertex but
4510167291SKenneth E. Jansen          !acos(0.75) is a pretty big angle. On coarse meshes this number may
4610167291SKenneth E. Jansen          !have to be dialed down to get normal on curved walls correct.
4710167291SKenneth E. Jansen
4810167291SKenneth E. Jansen          allocate(wNorm(nshg,nsd))
4910167291SKenneth E. Jansen          wNorm = 0
5010167291SKenneth E. Jansen
5110167291SKenneth E. Jansen          do iblk=1,nelblb !loop over all the bounday element blocks at current process
5210167291SKenneth E. Jansen            npro=lcblkb(1,iblk+1)-lcblkb(1,iblk)
5310167291SKenneth E. Jansen            nenbl=lcblkb(6,iblk)
5410167291SKenneth E. Jansen            nshl=lcblkb(9,iblk)
5510167291SKenneth E. Jansen            allocate(ienb(nshl))
5610167291SKenneth E. Jansen            do i=1,npro  ! loop over boundary elements
5710167291SKenneth E. Jansen!             isfID=miBCB(iblk)%p(i,2) ! miBCB(2) is the surfID
5810167291SKenneth E. Jansen              ienb(1:nshl)=mienb(iblk)%p(i,1:nshl)
5910167291SKenneth E. Jansen              do j=1,nenbl ! loop over boundary nodes
6010167291SKenneth E. Jansen                nn=ienb(j) !nn is the global index of suction node
6110167291SKenneth E. Jansen                if(j.ne.1.and.j.ne.nenbl)then
6210167291SKenneth E. Jansen                  aedg(:)=x(ienb(j+1),:)-x(nn,:)
6310167291SKenneth E. Jansen                  bedg(:)=x(ienb(j-1),:)-x(nn,:)
6410167291SKenneth E. Jansen                elseif(j.eq.1)then
6510167291SKenneth E. Jansen                  aedg(:)=x(ienb(j+1),:)-x(nn,:)
6610167291SKenneth E. Jansen                  bedg(:)=x(ienb(nenbl),:)-x(nn,:)
6710167291SKenneth E. Jansen                elseif(j.eq.nenbl)then
6810167291SKenneth E. Jansen                  aedg(:)=x(ienb(1),:)-x(nn,:)
6910167291SKenneth E. Jansen                  bedg(:)=x(ienb(j-1),:)-x(nn,:)
7010167291SKenneth E. Jansen                endif
7110167291SKenneth E. Jansen                eNorm(1)=aedg(2)*bedg(3)-aedg(3)*bedg(2)
7210167291SKenneth E. Jansen                eNorm(2)=aedg(3)*bedg(1)-aedg(1)*bedg(3)
7310167291SKenneth E. Jansen                eNorm(3)=aedg(1)*bedg(2)-aedg(2)*bedg(1)
7410167291SKenneth E. Jansen
7510167291SKenneth E. Jansen                call addElemNormal(wNorm, eNorm, nn)
7610167291SKenneth E. Jansen
7710167291SKenneth E. Jansen              enddo  !end do j = 1, nenbl
7810167291SKenneth E. Jansen            enddo  !end do i = 1, npro
7910167291SKenneth E. Jansen            deallocate(ienb)
8010167291SKenneth E. Jansen          enddo  !end do iblk = 1, nelblk
8110167291SKenneth E. Jansen
8210167291SKenneth E. Jansen          !Now communicate between parts. Note: there still may be a
8310167291SKenneth E. Jansen          !conflict if the boundary edge is on a corner, so a separate array
8410167291SKenneth E. Jansen          !needs to be used to allow us to repeat the z-norm comparison.
8510167291SKenneth E. Jansen          if(numpe.gt.1) then
8610167291SKenneth E. Jansen            allocate(wNormTmp(nshg,nsd))
8710167291SKenneth E. Jansen            wNormTmp = wNorm
8810167291SKenneth E. Jansen            call commu(wNormTmp(:,:), ilwork, 3, 'in ')
8910167291SKenneth E. Jansen            wNormTmp = wNormTmp - wNorm   !undo the addition in commu to just
9010167291SKenneth E. Jansen                                          !get contributions from other processors
9110167291SKenneth E. Jansen
9210167291SKenneth E. Jansen            do i = 1,nshg                 !loop over each node and add the off
9310167291SKenneth E. Jansen              wn(:) = wNormTmp(i,:)       !processor contribution. If this is
9410167291SKenneth E. Jansen              call addElemNormal(wNorm, wn, i) !an interior node, zero gets
9510167291SKenneth E. Jansen            enddo                              !added to zero.
9610167291SKenneth E. Jansen          endif
9710167291SKenneth E. Jansen
9810167291SKenneth E. Jansen          !It is unclear how to treat periodic boundary conditions, i.e.
9910167291SKenneth E. Jansen          !should periodic walls have zero norm, a correct normal for each
10010167291SKenneth E. Jansen          !side, or a correct normal for the master. With bc3per commented,
10110167291SKenneth E. Jansen          !the wall normal should be calculated as though there is no
10210167291SKenneth E. Jansen          !periodic BC.
10310167291SKenneth E. Jansen!         call bc3per(iBC,wNorm(:,:),iper,ilwork,3)
10410167291SKenneth E. Jansen
10510167291SKenneth E. Jansen          if(numpe.gt.1) call commu(wNorm(:,:),ilwork,3,'out')
10610167291SKenneth E. Jansen
10710167291SKenneth E. Jansen          do nn = 1,nshg
10810167291SKenneth E. Jansen            wn(:) = wNorm(nn,:)
10910167291SKenneth E. Jansen            mag = (wn(1)**2+wn(2)**2+wn(3)**2)**0.5
11010167291SKenneth E. Jansen
11110167291SKenneth E. Jansen            if(mag > 0) wNorm(nn,:) = wn(:)/mag
11210167291SKenneth E. Jansen          enddo
11310167291SKenneth E. Jansen
11410167291SKenneth E. Jansen          zwall=0.0762d00
11510167291SKenneth E. Jansen          do nn=1,nshg ! normalize wNorm
11610167291SKenneth E. Jansen            wn(:)=wNorm(nn,:)
11710167291SKenneth E. Jansen            mag=(wn(1)**2+wn(2)**2+wn(3)**2)**0.5
11810167291SKenneth E. Jansen            if(mag.ne.0) then
11910167291SKenneth E. Jansen              if(((dabs(dabs(x(nn,3))-zwall)).lt.0.000000001d00)
12010167291SKenneth E. Jansen     &                  .and. (dabs(wn(3)/mag) .lt.0.99)) then
12110167291SKenneth E. Jansen! logic to trap wall normals in z above failed, possibley in commu? so
12210167291SKenneth E. Jansen! report it and "fix it"
12310167291SKenneth E. Jansen                write(*,*)
12410167291SKenneth E. Jansen     &             "Warning: still screwing up wall Norm, proc ",
12510167291SKenneth E. Jansen     &              myrank, ", node ", nn,
12610167291SKenneth E. Jansen     &              ", x = (", x(nn,1), x(nn,2), x(nn,3), ")",
12710167291SKenneth E. Jansen     &              ", n = (", wn(1), wn(2), wn(3), ")"
12810167291SKenneth E. Jansen!                write(myrank+1000,1234)
12910167291SKenneth E. Jansen!     &              nn,x(nn,1),x(nn,2),x(nn,3),wn(1),wn(2),wn(3)
13010167291SKenneth E. Jansen                wn(1)=0
13110167291SKenneth E. Jansen                wn(2)=0
13210167291SKenneth E. Jansen                wn(3)=x(nn,3)
13310167291SKenneth E. Jansen                mag=abs(wn(3))
13410167291SKenneth E. Jansen              endif
13510167291SKenneth E. Jansen              wNorm(nn,:)=wn(:)/mag
13610167291SKenneth E. Jansen            endif
13710167291SKenneth E. Jansen          enddo
13810167291SKenneth E. Jansen
13910167291SKenneth E. Jansen1234      format(i5,6(2x,e12.5))
14010167291SKenneth E. Jansen!!DEBUG
14110167291SKenneth E. Jansen!!      dimension   u(nshg,n)
14210167291SKenneth E. Jansen!      do j = 1,3
14310167291SKenneth E. Jansen!        do i = 1,nshg
14410167291SKenneth E. Jansen!         u(i,j)=9.876543e21
14510167291SKenneth E. Jansen!          wallBCout(i, j  ) = BC(   i, j+2)
14610167291SKenneth E. Jansen!          wallBCout(i, j+3) = wNorm(i, j  )
14710167291SKenneth E. Jansen!        enddo
14810167291SKenneth E. Jansen!      enddo
14910167291SKenneth E. Jansen!      call write_restart(myrank,100002,nshg,6,wallBCout,wallBCout)
15010167291SKenneth E. Jansen!!DEBUG
15110167291SKenneth E. Jansen
15210167291SKenneth E. Jansen          wNormCalced = .true.
15310167291SKenneth E. Jansen
15410167291SKenneth E. Jansen          return
15513c4f35aSKenneth E. Jansen        end subroutine
15610167291SKenneth E. Jansen
15710167291SKenneth E. Jansen        subroutine addElemNormal(wNorm, eNorm, nn)
15810167291SKenneth E. Jansen        !Adds contributions from the wall element normal en into the wall
15910167291SKenneth E. Jansen        !normal array wNorm at index nn.
16010167291SKenneth E. Jansen
16110167291SKenneth E. Jansen          include "common.h"
16210167291SKenneth E. Jansen
16310167291SKenneth E. Jansen          real*8 :: wNorm(nshg, nsd)
16410167291SKenneth E. Jansen          real*8 :: eNorm(nsd), en(nsd), wn(nsd)
16510167291SKenneth E. Jansen          real*8 :: mag, oDotn
16610167291SKenneth E. Jansen          integer :: nn
16710167291SKenneth E. Jansen
16810167291SKenneth E. Jansen          wn(:)=wNorm(nn,:)
16910167291SKenneth E. Jansen          mag=(wn(1)**2+wn(2)**2+wn(3)**2)**0.5
17010167291SKenneth E. Jansen
17110167291SKenneth E. Jansen          if(mag < 1.0e-14) then  ! no direction yet so give it this
17210167291SKenneth E. Jansen            wNorm(nn,:) = eNorm(:)
17310167291SKenneth E. Jansen
17410167291SKenneth E. Jansen          else ! wNorm already has a direction
17510167291SKenneth E. Jansen            wn(:)=wn(:)/mag
17610167291SKenneth E. Jansen
17710167291SKenneth E. Jansen            mag=(eNorm(1)**2+eNorm(2)**2+eNorm(3)**2)**0.5
17810167291SKenneth E. Jansen            en(:)=eNorm(:)/mag
17910167291SKenneth E. Jansen
18010167291SKenneth E. Jansen            ! when we are in corners we need to establish a precidence since
18110167291SKenneth E. Jansen            ! averaging will produce irregular results.  In this first pass,
18210167291SKenneth E. Jansen            ! let's let the vector with the largest z-normal take precedence.
18310167291SKenneth E. Jansen            oDotn=en(1)*wn(1)+en(2)*wn(2)+en(3)*wn(3)
18410167291SKenneth E. Jansen            if(abs(oDotn) < cornerCosThresh) then  ! this is a corner-Don't average
18510167291SKenneth E. Jansen              if(abs(en(3)) > abs(wn(3))) then
18610167291SKenneth E. Jansen                 wNorm(nn,:) = eNorm(:)  ! overwrite wNorm
18710167291SKenneth E. Jansen              endif  ! since no else, wNorm unaffected by weaker z
18810167291SKenneth E. Jansen                     ! e.g., for the else, eNorm is ignored
18910167291SKenneth E. Jansen            else !  not a corner so average
19010167291SKenneth E. Jansen              wNorm(nn,:)=wNorm(nn,:) + eNorm(:)
19110167291SKenneth E. Jansen            endif  !end if(abs(oDotn) < cornerCosThresh)
19210167291SKenneth E. Jansen          endif  !end mag < 1e-14
19310167291SKenneth E. Jansen
19410167291SKenneth E. Jansen        end subroutine
19510167291SKenneth E. Jansen
19610167291SKenneth E. Jansen      end module
19710167291SKenneth E. Jansen
198