xref: /phasta/phSolver/common/genbc.f (revision 513954ef803c300cddba2bb96b4a5dac0b93489a)
159599516SKenneth E. Jansen      subroutine genBC (iBC,  BC,   x,   ilwork, iper)
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc----------------------------------------------------------------------
459599516SKenneth E. Jansenc  This routine generates the essential prescribed boundary conditions.
559599516SKenneth E. Jansenc
659599516SKenneth E. Jansenc input:
759599516SKenneth E. Jansenc  iBC   (nshg)        : boundary condition code
859599516SKenneth E. Jansenc  nBC   (nshg)        : boundary condition mapping array
959599516SKenneth E. Jansenc
1059599516SKenneth E. Jansenc output:
1159599516SKenneth E. Jansenc  BC    (nshg,ndofBC) : The constraint data for prescribed BC
1259599516SKenneth E. Jansenc
1359599516SKenneth E. Jansenc
1459599516SKenneth E. Jansenc Note: genBC1 reduces the input data for the velocity. In the
1559599516SKenneth E. Jansenc       case of varying velocity direction in the generation, the
1659599516SKenneth E. Jansenc       results may not be correct. (since a linearity assumption is
1759599516SKenneth E. Jansenc       made in the generation).
1859599516SKenneth E. Jansenc
1959599516SKenneth E. Jansenc
2059599516SKenneth E. Jansenc Farzin Shakib, Spring 1986.
2159599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
2259599516SKenneth E. Jansenc----------------------------------------------------------------------
2359599516SKenneth E. Jansenc
2459599516SKenneth E. Jansen      use slpw
2559599516SKenneth E. Jansen      use readarrays            ! used to access BCinp, nBC
2659599516SKenneth E. Jansen      use specialBC ! filling acs here
2759599516SKenneth E. Jansen      include "common.h"
2859599516SKenneth E. Jansenc
2959599516SKenneth E. Jansen      dimension iBC(nshg),                nsurf(nshg),
3059599516SKenneth E. Jansen     &            BC(nshg,ndofBC),
3159599516SKenneth E. Jansen     &            x(numnp,nsd),             ilwork(nlwork),
3259599516SKenneth E. Jansen     &            iper(nshg)
3359599516SKenneth E. Jansenc
3459599516SKenneth E. Jansenc BCinp for each point has:
3559599516SKenneth E. Jansenc   D T P c11 c12 c13 M1 c21 c22 c23 M2 theta S1 S2 S3...
3659599516SKenneth E. Jansenc   1 2 3 4   5   6   7  8   9   10  11 12    13 14 15...
3759599516SKenneth E. Jansenc Remember, ndof=nsd+2+nsclr
3859599516SKenneth E. Jansenc
3959599516SKenneth E. Jansenc Arrays in the following 1 line are now dimensioned in readnblk
4059599516SKenneth E. Jansenc        dimension BCinp(numpbc,ndof+7)
4159599516SKenneth E. Jansenc
4259599516SKenneth E. Jansen      dimension BCtmp(nshg,ndof+7)
4359599516SKenneth E. Jansenc
4459599516SKenneth E. Jansenc ndof+7= 3(thermos) + (nsd-1)*(nsd+1) + nscalars + 1 (theta)
4559599516SKenneth E. Jansenc                       #vect *(vec dir +mag)
4659599516SKenneth E. Jansenc
4759599516SKenneth E. Jansenc.... --------------------------->  Input  <---------------------------
4859599516SKenneth E. Jansenc
4959599516SKenneth E. Jansenc.... convert boundary condition data
5059599516SKenneth E. Jansenc
5159599516SKenneth E. Jansen      BCtmp = zero
5259599516SKenneth E. Jansenc
5359599516SKenneth E. Jansen      if(numpbc.ne.0) then
5459599516SKenneth E. Jansen         do i = 1, ndof+7
5559599516SKenneth E. Jansen            where (nBC(:) .ne. 0) BCtmp(:,i) = BCinp(nBC(:),i)
5659599516SKenneth E. Jansen         enddo
5759599516SKenneth E. Jansen         deallocate(BCinp)
5859599516SKenneth E. Jansen      endif
5959599516SKenneth E. Jansen
6059599516SKenneth E. Jansenc
6159599516SKenneth E. Jansen      if(any(BCtmp(:,12).ne.0)) then
6259599516SKenneth E. Jansen         iabc=1
6359599516SKenneth E. Jansen         allocate (acs(nshg,2))
6459599516SKenneth E. Jansen         where (btest(iBC,10))
6559599516SKenneth E. Jansen            acs(:,1) = cos(BCtmp(:,12))
6659599516SKenneth E. Jansen            acs(:,2) = sin(BCtmp(:,12))
6759599516SKenneth E. Jansen         endwhere
6859599516SKenneth E. Jansen      endif
6959599516SKenneth E. Jansen
7059599516SKenneth E. Jansenc
7159599516SKenneth E. Jansenc.... ----------------------> Wall Normals  <--------------------------
7259599516SKenneth E. Jansenc (calculate the normal and adjust BCinp to the true normal as needed)
7359599516SKenneth E. Jansenc
7459599516SKenneth E. Jansen      if(navier.eq.1)then
7559599516SKenneth E. Jansen         call genwnm (iBC,  BCtmp,   x,   ilwork, iper, nsurf)
7659599516SKenneth E. Jansen      endif
7759599516SKenneth E. Jansenc  determine the first point off the wall for each wall node
7859599516SKenneth E. Jansen      if(navier.eq.1)then
7959599516SKenneth E. Jansen         call genotwn (x,BCtmp, iBC, nsurf)
8059599516SKenneth E. Jansen      endif
8159599516SKenneth E. Jansenc.... ------------------------>  Conversion  <-------------------------
8259599516SKenneth E. Jansenc
8359599516SKenneth E. Jansenc.... convert the input boundary conditions to condensed version
8459599516SKenneth E. Jansenc
8559599516SKenneth E. Jansen      BC = zero
8659599516SKenneth E. Jansenc
87*513954efSKenneth E. Jansen      if(myrank.eq.0) write(*,*) 'Navier is set to ', navier
8859599516SKenneth E. Jansen      if(navier.eq.1)then ! zero navier means Euler simulation
8959599516SKenneth E. Jansen         call genBC1 (BCtmp,  iBC,  BC)
90*513954efSKenneth E. Jansen      elseif(matflg(1,1).eq.0)then !  compressible code
9159599516SKenneth E. Jansen         allocate(BCtmpg(nshg,ndof+7))
9259599516SKenneth E. Jansen         allocate(BCg(nshg,ndofBC))
9359599516SKenneth E. Jansen         allocate(iBCg(nshg))
9459599516SKenneth E. Jansen         BCtmpg=BCtmp
9559599516SKenneth E. Jansen         iBCg=iBC
9659599516SKenneth E. Jansen         call genBC1 (BCtmp,  iBC,  BC)
9759599516SKenneth E. Jansenc... genBC1 convert BCtmp to BC
9859599516SKenneth E. Jansen         BCg=BC
9959599516SKenneth E. Jansen         icdg=icd
10059599516SKenneth E. Jansenc... find slip wall
10159599516SKenneth E. Jansen         call findslpw(x,ilwork,iper,iBC)
10259599516SKenneth E. Jansenc... apply slip wall condition to wall nodes
10359599516SKenneth E. Jansen         do i=1,nslpwnd
10459599516SKenneth E. Jansen            nn=idxslpw(i)
10559599516SKenneth E. Jansen            call slpwBC(mpslpw(nn,1),mpslpw(nn,2),iBCg(nn),
10659599516SKenneth E. Jansen     &               BCg(nn,:),  BCtmpg(nn,:),
10759599516SKenneth E. Jansen     &               iBC(nn),    BC(nn,:),
10859599516SKenneth E. Jansen     &               wlnorm(nn,:,:)                     )
10959599516SKenneth E. Jansen         enddo
11059599516SKenneth E. Jansen         icd=icdg
111*513954efSKenneth E. Jansen      else
112*513954efSKenneth E. Jansen         if(myrank.eq.0) write(*,*) 'Incompressible code not able to do inviscid at this time'
11359599516SKenneth E. Jansen      endif
11459599516SKenneth E. Jansen
11559599516SKenneth E. Jansen
11659599516SKenneth E. Jansenc
11759599516SKenneth E. Jansenc.... --------------------------->  Echo  <----------------------------
11859599516SKenneth E. Jansenc
11959599516SKenneth E. Jansenc.... echo the input data
12059599516SKenneth E. Jansenc
12159599516SKenneth E. Jansen      if (necho .lt. 3) then
12259599516SKenneth E. Jansen         nn = 0
12359599516SKenneth E. Jansen         do n = 1, nshg
12459599516SKenneth E. Jansen            if (nBC(n) .ne. 0) then
12559599516SKenneth E. Jansen               nn = nn + 1
12659599516SKenneth E. Jansen               if(mod(nn,50).eq.1)
12759599516SKenneth E. Jansen     &              write(iecho,1000)ititle,(j,j=1,ndofBC)
12859599516SKenneth E. Jansen               write (iecho,1100) n, (BC(n,i),i=1,ndofBC)
12959599516SKenneth E. Jansen            endif
13059599516SKenneth E. Jansen         enddo
13159599516SKenneth E. Jansen      endif
13259599516SKenneth E. Jansenc
13359599516SKenneth E. Jansenc.... return
13459599516SKenneth E. Jansenc
13559599516SKenneth E. Jansen      return
13659599516SKenneth E. Jansenc
13759599516SKenneth E. Jansen 1000 format(a80,//,
13859599516SKenneth E. Jansen     &' P r e s c r i b e d   B o u n d a r y   C o n d i t i o n s',//,
13959599516SKenneth E. Jansen     &  '    Node  ',/,
14059599516SKenneth E. Jansen     &  '   Number ',5x,6('BC',i1,:,10x))
14159599516SKenneth E. Jansen 1100 format(1p,2x,i5,3x,6(e12.5,1x))
14259599516SKenneth E. Jansenc
14359599516SKenneth E. Jansen      end
14459599516SKenneth E. Jansen
14559599516SKenneth E. Jansen
14659599516SKenneth E. Jansen
14759599516SKenneth E. Jansen
14859599516SKenneth E. Jansen
14959599516SKenneth E. Jansen
15059599516SKenneth E. Jansen
15159599516SKenneth E. Jansen
15259599516SKenneth E. Jansen
15359599516SKenneth E. Jansen
15459599516SKenneth E. Jansen      subroutine genwnm (iBC, BCtmp,    x,   ilwork, iper, nsurf)
15559599516SKenneth E. Jansenc----------------------------------------------------------------------
15659599516SKenneth E. Jansenc  This routine generates the normal to a wall
15759599516SKenneth E. Jansenc
15859599516SKenneth E. Jansenc input:
15959599516SKenneth E. Jansenc  iBC   (nshg)        : boundary condition code
16059599516SKenneth E. Jansenc  nBC   (nshg)        : boundary condition mapping array
16159599516SKenneth E. Jansenc
16259599516SKenneth E. Jansenc output:
16359599516SKenneth E. Jansenc  BCtmp    (nshg,ndof+6) : The constraint data for prescribed BC
16459599516SKenneth E. Jansenc
16559599516SKenneth E. Jansenc
16659599516SKenneth E. Jansenc----------------------------------------------------------------------
16759599516SKenneth E. Jansenc
16859599516SKenneth E. Jansen      use turbSA
16959599516SKenneth E. Jansen      use pointer_data          ! used for mienb, mibcb
17059599516SKenneth E. Jansen      include "common.h"
17159599516SKenneth E. Jansen      include "mpif.h"
17259599516SKenneth E. Jansenc
17359599516SKenneth E. Jansen      character*20 fname1,  fmt1
17459599516SKenneth E. Jansen      character*5  cname
17559599516SKenneth E. Jansen      dimension iBC(nshg),                  iper(nshg),
17659599516SKenneth E. Jansen     &            x(numnp,nsd),             ilwork(nlwork)
17759599516SKenneth E. Jansenc
17859599516SKenneth E. Jansen      dimension  BCtmpSAV(nshg,ndof+7)
17959599516SKenneth E. Jansen      dimension  BCtmp(nshg,ndof+7),      fBC(nshg,ndofBC),
18059599516SKenneth E. Jansen     &            e1(3),                    e2(3),
18159599516SKenneth E. Jansen     &            elnrm(3),                 asum(numnp)
18259599516SKenneth E. Jansenc
18359599516SKenneth E. Jansen      integer sid, nsidg
18459599516SKenneth E. Jansen      integer nsurf(nshg), ivec(nsd)
18559599516SKenneth E. Jansen      logical :: firstvisit(nshg)
18659599516SKenneth E. Jansen      real*8 BCvecs(2,nsd)
18759599516SKenneth E. Jansen      integer, allocatable :: ienb(:)
18859599516SKenneth E. Jansen      dimension wnmdb(nshg,nsd)
18959599516SKenneth E. Jansenc
19059599516SKenneth E. Jansenc  wnrm is dimensioned nshg but the math is only done for straight sided
19159599516SKenneth E. Jansenc  elements at this point so wnrm will not be calculated for the hierarchic
19259599516SKenneth E. Jansenc  modes.  Note that the wall model creates a p.w. linear representation
19359599516SKenneth E. Jansenc  only at this time.
19459599516SKenneth E. Jansenc
19559599516SKenneth E. Jansen      allocate ( wnrm(nshg,3) )
19659599516SKenneth E. Jansenc
19759599516SKenneth E. Jansenc.... ----------------------> Wall Normals  <--------------------------
19859599516SKenneth E. Jansenc (calculate the normal and adjust BCinp to the true normal as needed)
19959599516SKenneth E. Jansenc
20059599516SKenneth E. Jansenc
20159599516SKenneth E. Jansen      asum = zero
20259599516SKenneth E. Jansen      wnrm = zero
20359599516SKenneth E. Jansen
20459599516SKenneth E. Jansenc
20559599516SKenneth E. Jansenc....  Save a copy of BCtmp so that after we calculate the normals we
20659599516SKenneth E. Jansenc      can recover the comp3 information.
20759599516SKenneth E. Jansenc
20859599516SKenneth E. Jansen      BCtmpSAV=BCtmp
20959599516SKenneth E. Jansenc
21059599516SKenneth E. Jansenc Count out the number of surface ID's on-processor, and map them
21159599516SKenneth E. Jansen      call gensidcount(nsidg)
21259599516SKenneth E. Jansen
21359599516SKenneth E. Jansenc
21459599516SKenneth E. Jansen      if(nsidg.gt.0) then       ! if there are any surfID's
21559599516SKenneth E. Jansen         nsurf(:) = 0
21659599516SKenneth E. Jansen         do k = 1, nsidg        ! loop over Surface ID's
21759599516SKenneth E. Jansen            sid = sidmapg(k)
21859599516SKenneth E. Jansen            firstvisit(:)=.true.
21959599516SKenneth E. Jansen            wnrm(:,1:3)=zero
22059599516SKenneth E. Jansen            do iblk=1, nelblb   ! loop over boundary element blocks
22159599516SKenneth E. Jansen               npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
22259599516SKenneth E. Jansen               nenbl = lcblkb(6,iblk)
22359599516SKenneth E. Jansen               nshl = lcblkb(9,iblk)
22459599516SKenneth E. Jansen               allocate( ienb(nshl) )
22559599516SKenneth E. Jansen               do i = 1, npro   ! loop over boundary elements
22659599516SKenneth E. Jansen                  iBCB1=miBCB(iblk)%p(i,1)
22759599516SKenneth E. Jansen                  iBCB2=miBCB(iblk)%p(i,2)
22859599516SKenneth E. Jansen                  ienb(1:nshl)=mienb(iblk)%p(i,1:nshl)
22959599516SKenneth E. Jansenc don't bother with elements that aren't on modeled surfaces
23059599516SKenneth E. Jansenc              if ( not          (wall set  and   traction set)    )
23159599516SKenneth E. Jansen                  if (.not.(btest(iBCB1,2).and.btest(iBCB1,4)))
23259599516SKenneth E. Jansen     &                 cycle
23359599516SKenneth E. Jansenc don't bother with elements that don't lie on the current surface
23459599516SKenneth E. Jansen                  if (iBCB2.ne.sid) cycle
23559599516SKenneth E. Jansenc
23659599516SKenneth E. Jansenc.... calculate this element's area-weighted normal vector
23759599516SKenneth E. Jansenc
23859599516SKenneth E. Jansen                  e1 = x(ienb(2),:)-x(ienb(1),:)
23959599516SKenneth E. Jansen                  e2 = x(ienb(3),:)-x(ienb(1),:)
24059599516SKenneth E. Jansen                  elnrm(1) = e1(2)*e2(3)-e1(3)*e2(2)
24159599516SKenneth E. Jansen                  elnrm(2) = e1(3)*e2(1)-e1(1)*e2(3)
24259599516SKenneth E. Jansen                  elnrm(3) = e1(1)*e2(2)-e1(2)*e2(1)
24359599516SKenneth E. Jansenc Tetrahedral elements have negative volumes in phastaI, so
24459599516SKenneth E. Jansenc the normals calculated from the boundary faces must be inverted
24559599516SKenneth E. Jansenc to point into the fluid
24659599516SKenneth E. Jansen                  if(nenbl.eq.3) elnrm(:)=-elnrm(:)
24759599516SKenneth E. Jansenc
24859599516SKenneth E. Jansenc.... add area-weighted normals to the nodal tallies for this surface
24959599516SKenneth E. Jansenc
25059599516SKenneth E. Jansen                  do j = 1, nenbl ! loop over elt boundary nodes
25159599516SKenneth E. Jansen                     nn=ienb(j) ! global node number
25259599516SKenneth E. Jansen                     if(firstvisit(nn)) then
25359599516SKenneth E. Jansen                        firstvisit(nn)=.false.
25459599516SKenneth E. Jansen                        nsurf(nn)=nsurf(nn)+1
25559599516SKenneth E. Jansen                        if(nsurf(nn).eq.1) BCtmp(nn,4:6)=zero
25659599516SKenneth E. Jansen                        if(nsurf(nn).eq.2) BCtmp(nn,8:10)=zero
25759599516SKenneth E. Jansen                     endif
25859599516SKenneth E. Jansen                     wnrm(nn,:)=wnrm(nn,:)+elnrm
25959599516SKenneth E. Jansen                  enddo         ! loop over elt boundary nodes
26059599516SKenneth E. Jansen               enddo            ! end loop over boundary elements in block
26159599516SKenneth E. Jansen               deallocate(ienb)
26259599516SKenneth E. Jansen            enddo               ! end loop over boundary element blocks
26359599516SKenneth E. Jansenc Now we have all of this surface's contributions to wall normals
26459599516SKenneth E. Jansenc for all nodes, along with an indication of how many surfaces
26559599516SKenneth E. Jansenc each node has encountered so far.  Contributions from other processors
26659599516SKenneth E. Jansenc should now be accumulated for this surface
26759599516SKenneth E. Jansenc
26859599516SKenneth E. Jansenc axisymm BC's need BC (and we have not built it yet) so we need to create
26959599516SKenneth E. Jansenc the entries it needs.
27059599516SKenneth E. Jansenc
27159599516SKenneth E. Jansen            if(iabc==1) then
27259599516SKenneth E. Jansen               where (btest(iBC,10))
27359599516SKenneth E. Jansen                  fBC(:,1) = cos(BCtmp(:,12))
27459599516SKenneth E. Jansen                  fBC(:,2) = sin(BCtmp(:,12))
27559599516SKenneth E. Jansen               endwhere
27659599516SKenneth E. Jansen
27759599516SKenneth E. Jansenc before the commu we need to rotate the residual vector for axisymmetric
27859599516SKenneth E. Jansenc boundary conditions (so that off processor periodicity is a dof add instead
27959599516SKenneth E. Jansenc of a dof combination).  Take care of all nodes now so periodicity, like
28059599516SKenneth E. Jansenc commu is a simple dof add.
28159599516SKenneth E. Jansen               call rotabc(wnrm, iBC, 'in ')
28259599516SKenneth E. Jansen            endif
28359599516SKenneth E. Jansen
28459599516SKenneth E. Jansen            if (numpe.gt.1) then
28559599516SKenneth E. Jansenc.... add areas and norms contributed from other processors
28659599516SKenneth E. Jansen               call commu (wnrm, ilwork, 3, 'in ')
28759599516SKenneth E. Jansen            endif
28859599516SKenneth E. Jansenc.... account for periodicity and axisymmetry
28959599516SKenneth E. Jansen            call bc3per(iBC,wnrm, iper,ilwork, 3)
29059599516SKenneth E. Jansenc at this point the master has all the information (slaves are zeroed and
29159599516SKenneth E. Jansenc waiting to be told what the master has...lets do that).
29259599516SKenneth E. Jansenc
29359599516SKenneth E. Jansenc.... local periodic (and axisymmetric) boundary conditions (no communications)
29459599516SKenneth E. Jansen            wnmdb=wnrm
29559599516SKenneth E. Jansen            do i = 1,numnp      ! only use the vertices in the normal calc
29659599516SKenneth E. Jansen               wnrm(i,:) = wnrm(iper(i),:)
29759599516SKenneth E. Jansen            enddo
29859599516SKenneth E. Jansen            wnmdb=wnrm
29959599516SKenneth E. Jansen            if (numpe.gt.1) then
30059599516SKenneth E. Jansenc.... tell other processors the resulting (and now complete) sums
30159599516SKenneth E. Jansen               call commu (wnrm, ilwork, 3, 'out')
30259599516SKenneth E. Jansen            endif
30359599516SKenneth E. Jansenc Axisymmetric slaves need to be rotated back
30459599516SKenneth E. Jansen            if(iabc==1) then    !are there any axisym bc's
30559599516SKenneth E. Jansen               call rotabc(wnrm, iBC, 'out')
30659599516SKenneth E. Jansen            endif
30759599516SKenneth E. Jansenc Now all nodes have all the normal contributions for this surface,
30859599516SKenneth E. Jansenc including those from off-processor and periodic nodes.
30959599516SKenneth E. Jansenc We distribute these summed vectors to the proper slots in BCtmp,
31059599516SKenneth E. Jansenc according to how many surfaces each node has seen so far
31159599516SKenneth E. Jansen            do nn = 1, nshg
31259599516SKenneth E. Jansen               if(nsurf(nn).eq.1)
31359599516SKenneth E. Jansen     &              BCtmp(nn,4:6)=BCtmp(nn,4:6)+wnrm(nn,:)
31459599516SKenneth E. Jansen               if(nsurf(nn).eq.2)
31559599516SKenneth E. Jansen     &              BCtmp(nn,8:10)=BCtmp(nn,8:10)+wnrm(nn,:)
31659599516SKenneth E. Jansen               if(nsurf(nn).gt.2)
31759599516SKenneth E. Jansen     &              BCtmp(nn,4:6)=BCtmp(nn,4:6)+wnrm(nn,:)
31859599516SKenneth E. Jansen            enddo
31959599516SKenneth E. Jansenc That's all for this surface; move on to the next
32059599516SKenneth E. Jansen         enddo                  ! end loop over surface ID's
3216b966dd8SCameron Smith         deallocate( sidmapg )
32259599516SKenneth E. Jansen
32359599516SKenneth E. Jansenc
32459599516SKenneth E. Jansenc.... complete the averaging, adjust iBC, adjust BCtmp
32559599516SKenneth E. Jansenc
32659599516SKenneth E. Jansen         wnrm(:,1:3)=zero
32759599516SKenneth E. Jansen         do nn = 1, numnp       ! loop over all nodes
32859599516SKenneth E. Jansenc leave nodes without wall-modeled surfaces unchanged
32959599516SKenneth E. Jansen            if(nsurf(nn).eq.0) cycle
33059599516SKenneth E. Jansenc
33159599516SKenneth E. Jansenc mark this as a node with comp3
33259599516SKenneth E. Jansenc
33359599516SKenneth E. Jansen            ikp=0
33459599516SKenneth E. Jansen            if(ibits(iBC(nn),3,3).eq.7 ) ikp=1
33559599516SKenneth E. Jansenc         if(ibits(ibc(nn),3,3).eq.7 .and. BCtmp(nn,7).eq.zero) cycle
33659599516SKenneth E. Jansenc find out which velocity BC's were set by the user, and cancel them
33759599516SKenneth E. Jansen            ixset=ibits(iBC(nn),3,1)
33859599516SKenneth E. Jansen            iyset=ibits(iBC(nn),4,1)
33959599516SKenneth E. Jansen            izset=ibits(iBC(nn),5,1)
34059599516SKenneth E. Jansen            iBC(nn)=iBC(nn)-ixset*8-iyset*16-izset*32
34159599516SKenneth E. Jansenc
34259599516SKenneth E. Jansenc save this stripped iBC for later un-fixing
34359599516SKenneth E. Jansenc
34459599516SKenneth E. Jansen            iBCSAV=iBC(nn)
34559599516SKenneth E. Jansenc
34659599516SKenneth E. Jansen            if(abs(itwmod).eq.1) then ! slip velocity wall-model
34759599516SKenneth E. Jansenc If we're using the slip-velocity wall-model, then the velocity
34859599516SKenneth E. Jansenc boundary condition for this node will depend on how many modeled
34959599516SKenneth E. Jansenc surfaces share it...
35059599516SKenneth E. Jansen               if(nsurf(nn).eq.1) then ! 1 modeled wall
35159599516SKenneth E. Jansenc   If this node is part of only one modeled wall, then the component
35259599516SKenneth E. Jansenc   of the velocity normal to the surface is set to zero.  In this case
35359599516SKenneth E. Jansenc   only the first vector of BCtmp received normal contributions
35459599516SKenneth E. Jansen                  iBC(nn)=iBC(nn)+8 ! assume normal is x-dominated first
35559599516SKenneth E. Jansen                  wnrm(nn,:)=BCtmp(nn,4:6)
35659599516SKenneth E. Jansen                  if(abs(wnrm(nn,3)).gt.abs(wnrm(nn,2)))then ! z beats y
35759599516SKenneth E. Jansen                     if(abs(wnrm(nn,3)).gt.abs(wnrm(nn,1)))then ! z beats x
35859599516SKenneth E. Jansen                        iBC(nn)=iBC(nn)+24
35959599516SKenneth E. Jansen                     endif      ! z beats x
36059599516SKenneth E. Jansen                  endif         ! z beats y
36159599516SKenneth E. Jansen                  if(abs(wnrm(nn,2)).ge.abs(wnrm(nn,3)))then ! y beats z
36259599516SKenneth E. Jansen                     if(abs(wnrm(nn,2)).gt.abs(wnrm(nn,1)))then ! y beats x
36359599516SKenneth E. Jansen                        iBC(nn)=iBC(nn)+8
36459599516SKenneth E. Jansen                     endif      ! y beats x
36559599516SKenneth E. Jansen                  endif         ! y beats z           !(adjusted iBC)
36659599516SKenneth E. Jansen                  BCtmp(nn,7)=zero
36759599516SKenneth E. Jansen               else if(nsurf(nn).eq.2) then
36859599516SKenneth E. Jansenc   If this node is shared by two modeled walls, then each wall
36959599516SKenneth E. Jansenc   provides a normal vector along which the velocity must be zero.
37059599516SKenneth E. Jansenc   This leaves only one "free" direction, along which the flow is nonzero.
37159599516SKenneth E. Jansenc   The two normal vectors define a plane.  Any pair of non-parallel
37259599516SKenneth E. Jansenc   vectors in this plane can also be specified, leaving the same free
37359599516SKenneth E. Jansenc   direction.  Area-weighted average normal vectors for the two surfaces
37459599516SKenneth E. Jansenc   sharing this node have been stored in BCtmp(4:6) and BCtmp(8:10)
37559599516SKenneth E. Jansenc   We normalize the first of these, and then orthogonalize the second
37659599516SKenneth E. Jansenc   against the first.  After then normalizing the second, we choose which
37759599516SKenneth E. Jansenc   cartesian direction dominates each of the vectors, and adjust iBC
37859599516SKenneth E. Jansenc   and BCtmp accordingly
37959599516SKenneth E. Jansen                  e1=BCtmp(nn,4:6)
38059599516SKenneth E. Jansen                  wmag=e1(1)*e1(1)+e1(2)*e1(2)+e1(3)*e1(3)
38159599516SKenneth E. Jansen                  wmag=1./sqrt(wmag)
38259599516SKenneth E. Jansen                  e1=e1*wmag    ! first vector is normalized
38359599516SKenneth E. Jansenc
38459599516SKenneth E. Jansen                  e2=BCtmp(nn,8:10)
38559599516SKenneth E. Jansen                  wmag=e2(1)*e1(1)+e2(2)*e1(2)+e2(3)*e1(3) ! wmag=e2.n1
38659599516SKenneth E. Jansen                  e2=e2-wmag*e1 ! second vector is orthogonalized
38759599516SKenneth E. Jansenc
38859599516SKenneth E. Jansen                  wmag=e2(1)*e2(1)+e2(2)*e2(2)+e2(3)*e2(3)
38959599516SKenneth E. Jansen                  wmag=1./sqrt(wmag)
39059599516SKenneth E. Jansen                  e2=e2*wmag    ! second vector is normalized
39159599516SKenneth E. Jansenc
39259599516SKenneth E. Jansen               ! idom tells us which component is currently dominant
39359599516SKenneth E. Jansen               ! ivec(n) tells us which vector is dominated by comp n
39459599516SKenneth E. Jansen                  ivec(:)=0     ! initialize
39559599516SKenneth E. Jansen                  idom=1        ! assume x-comp dominates e1
39659599516SKenneth E. Jansen                  bigcomp=abs(e1(1))
39759599516SKenneth E. Jansen                  ivec(idom)=1
39859599516SKenneth E. Jansen                  do i = 2, nsd
39959599516SKenneth E. Jansen                     if(abs(e1(i)).gt.bigcomp) then
40059599516SKenneth E. Jansen                        ivec(idom)=0
40159599516SKenneth E. Jansen                        bigcomp=abs(e1(i))
40259599516SKenneth E. Jansen                        idom=i
40359599516SKenneth E. Jansen                        ivec(i)=1
40459599516SKenneth E. Jansen                     endif
40559599516SKenneth E. Jansen                  enddo
40659599516SKenneth E. Jansen                  if(idom.ne.1) then
40759599516SKenneth E. Jansen                     idom=1     ! assume x-comp dominates e2
40859599516SKenneth E. Jansen                  else
40959599516SKenneth E. Jansen                     idom=3     ! assume z-comp dominates e2
41059599516SKenneth E. Jansen                  endif
41159599516SKenneth E. Jansen                  bigcomp=abs(e2(idom))
41259599516SKenneth E. Jansen
41359599516SKenneth E. Jansen                  ivec(idom)=2
41459599516SKenneth E. Jansen                  do i = 1, nsd
41559599516SKenneth E. Jansen                     if(abs(e2(i)).gt.bigcomp) then
41659599516SKenneth E. Jansen                        if(ivec(i).eq.1) cycle
41759599516SKenneth E. Jansen                        ivec(idom)=0
41859599516SKenneth E. Jansen                        bigcomp=abs(e2(i))
41959599516SKenneth E. Jansen                        idom=i
42059599516SKenneth E. Jansen                        ivec(i)=2
42159599516SKenneth E. Jansen                     endif
42259599516SKenneth E. Jansen                  enddo
42359599516SKenneth E. Jansen               ! now we know which components dominate each vector
42459599516SKenneth E. Jansen                  ixset=0       !
42559599516SKenneth E. Jansen                  iyset=0       ! initialize
42659599516SKenneth E. Jansen                  izset=0       !
42759599516SKenneth E. Jansen                  if(ivec(1).ne.0) ixset=1
42859599516SKenneth E. Jansen                  if(ivec(2).ne.0) iyset=1
42959599516SKenneth E. Jansen                  if(ivec(3).ne.0) izset=1
43059599516SKenneth E. Jansen                  ncomp=ixset+iyset+izset
43159599516SKenneth E. Jansen                  if(ncomp.ne.2) write(*,*) 'WARNING: TROUBLE IN GENBC'
43259599516SKenneth E. Jansen                  BCvecs(1,1:3)=e1(:)
43359599516SKenneth E. Jansen                  BCvecs(2,1:3)=e2(:)
43459599516SKenneth E. Jansen                  if((ixset.eq.1).and.(iyset.eq.1)) then
43559599516SKenneth E. Jansen                     BCtmp(nn,4:6)=BCvecs(ivec(1),1:3)
43659599516SKenneth E. Jansen                     BCtmp(nn,8:10)=BCvecs(ivec(2),1:3)
43759599516SKenneth E. Jansen                  else if((ixset.eq.1).and.(izset.eq.1)) then
43859599516SKenneth E. Jansen                     BCtmp(nn,4:6)=BCvecs(ivec(1),1:3)
43959599516SKenneth E. Jansen                     BCtmp(nn,8:10)=BCvecs(ivec(3),1:3)
44059599516SKenneth E. Jansen                  else if((iyset.eq.1).and.(izset.eq.1)) then
44159599516SKenneth E. Jansen                     BCtmp(nn,4:6)=BCvecs(ivec(2),1:3)
44259599516SKenneth E. Jansen                     BCtmp(nn,8:10)=BCvecs(ivec(3),1:3)
44359599516SKenneth E. Jansen                  endif
44459599516SKenneth E. Jansen                  iBC(nn)=iBC(nn)+ixset*8+iyset*16+izset*32
44559599516SKenneth E. Jansen                  BCtmp(nn,7)=zero
44659599516SKenneth E. Jansen                  BCtmp(nn,11)=zero
44759599516SKenneth E. Jansen                  wnrm(nn,:)=BCtmp(nn,4:6)+BCtmp(nn,8:10)
44859599516SKenneth E. Jansen               else if(nsurf(nn).gt.2) then
44959599516SKenneth E. Jansenc     If this node is shared by more than two modeled walls, then
45059599516SKenneth E. Jansenc     it is a corner node and it will be treated as having no slip
45159599516SKenneth E. Jansenc     The normals for all surfaces beyond the first two were
45259599516SKenneth E. Jansenc     contributed to the first vector of BCtmp
45359599516SKenneth E. Jansen                  wnrm(nn,:)=BCtmp(nn,4:6)+BCtmp(nn,8:10)
45459599516SKenneth E. Jansen                  iBC(nn)=iBC(nn)+56
45559599516SKenneth E. Jansen                  BCtmp(nn,7)=zero
45659599516SKenneth E. Jansen               endif            ! curved surface
45759599516SKenneth E. Jansen            else if(abs(itwmod).eq.2) then ! effective viscosity wall-model
45859599516SKenneth E. Jansenc Otherwise, we're using the effective-viscosity wall-model and the
45959599516SKenneth E. Jansenc nodes on modeled surfaces are treated as no-slip
46059599516SKenneth E. Jansen               iBC(nn)=iBC(nn)+56 ! set 3-comp
46159599516SKenneth E. Jansen               if(itwmod.eq.-2) BCtmp(nn,7)=zero ! no slip
46259599516SKenneth E. Jansen               if(nsurf(nn).eq.1)
46359599516SKenneth E. Jansen     &              wnrm(nn,:)=BCtmp(nn,4:6)
46459599516SKenneth E. Jansen               if(nsurf(nn).ge.2)
46559599516SKenneth E. Jansen     &              wnrm(nn,:)=BCtmp(nn,4:6)+BCtmp(nn,8:10)
46659599516SKenneth E. Jansen            endif
46759599516SKenneth E. Jansenc Now normalize the wall normal for this node
46859599516SKenneth E. Jansen            if(itwmod.ne.0) then
46959599516SKenneth E. Jansen               wmag=sqrt(wnrm(nn,1)*wnrm(nn,1)
47059599516SKenneth E. Jansen     &              +wnrm(nn,2)*wnrm(nn,2)+wnrm(nn,3)*wnrm(nn,3))
47159599516SKenneth E. Jansen               wnrm(nn,:)=wnrm(nn,:)/wmag
47259599516SKenneth E. Jansen            endif
47359599516SKenneth E. Jansenc
47459599516SKenneth E. Jansenc.... put back the comp3 info for bctmp
47559599516SKenneth E. Jansenc
47659599516SKenneth E. Jansen            if(ikp.eq.1) then
47759599516SKenneth E. Jansen               iBC(nn)=iBCSAV+56 ! put it back to a comp3
47859599516SKenneth E. Jansen               BCtmp(nn,:)=BCtmpSAV(nn,:) ! ditto
47959599516SKenneth E. Jansen            endif
48059599516SKenneth E. Jansen         enddo                  ! loop over all nodes
48159599516SKenneth E. Jansen      endif                     ! end "if there are any surfID's"
48259599516SKenneth E. Jansenc
48359599516SKenneth E. Jansenc  If you are using the turbulence wall with axisymmetry you
48459599516SKenneth E. Jansenc  need to modify the axisymmetry angle to account for the discrete
48559599516SKenneth E. Jansenc  normals at the wall being different than the exact normals
48659599516SKenneth E. Jansenc
48759599516SKenneth E. Jansenc find the my normal, my partners normal and correct the angle
48859599516SKenneth E. Jansenc$$$
48959599516SKenneth E. Jansenc$$$        do i=1,numnp
49059599516SKenneth E. Jansenc$$$           wmag = wnrm(i,1) * wnrm(i,1)
49159599516SKenneth E. Jansenc$$$     &          + wnrm(i,2) * wnrm(i,2)
49259599516SKenneth E. Jansenc$$$     &          + wnrm(i,3) * wnrm(i,3)
49359599516SKenneth E. Jansenc$$$c
49459599516SKenneth E. Jansenc$$$c  only "fix" wall nodes...other nodes still have a zero in wnrm
49559599516SKenneth E. Jansenc$$$c
49659599516SKenneth E. Jansenc$$$           if((btest(iBC(i),12)).and.(wmag.ne.zero)) then
49759599516SKenneth E. Jansenc$$$              BCtmp(i,1)=acos( wnrm(i,1) * wnrm(iper(i),1)
49859599516SKenneth E. Jansenc$$$     &                       + wnrm(i,2) * wnrm(iper(i),2)
49959599516SKenneth E. Jansenc$$$     &                       + wnrm(i,3) * wnrm(iper(i),3) )
50059599516SKenneth E. Jansenc$$$           endif
50159599516SKenneth E. Jansenc$$$        enddo
50259599516SKenneth E. Jansenc
50359599516SKenneth E. Jansenc.... return
50459599516SKenneth E. Jansenc
50559599516SKenneth E. Jansen      return
50659599516SKenneth E. Jansenc
50759599516SKenneth E. Jansen      end
50859599516SKenneth E. Jansen
50959599516SKenneth E. Jansen
51059599516SKenneth E. Jansen      subroutine gensidcount(nsidg)
51159599516SKenneth E. Jansenc---------------------------------------------------------------------
51259599516SKenneth E. Jansenc
51359599516SKenneth E. Jansenc This routine counts up the total number of surface ID's across
51459599516SKenneth E. Jansenc all processors and makes a list of them
51559599516SKenneth E. Jansenc
51659599516SKenneth E. Jansenc Inputs:
51759599516SKenneth E. Jansenc     iBCB        natural boundary condition switches and surfIDs
51859599516SKenneth E. Jansenc
51959599516SKenneth E. Jansenc Outputs:
52059599516SKenneth E. Jansenc     nsidg       number of surface ID's globally (including all procs)
52159599516SKenneth E. Jansenc     sidmapg     global list of surface ID's, lowest to highest
52259599516SKenneth E. Jansenc
52359599516SKenneth E. Jansenc---------------------------------------------------------------------
52459599516SKenneth E. Jansen      use pointer_data ! access to miBCB
52559599516SKenneth E. Jansen      use turbSA ! access to sidmapg
52659599516SKenneth E. Jansen      include "common.h"
52759599516SKenneth E. Jansen      include "mpif.h"
52859599516SKenneth E. Jansenc
52959599516SKenneth E. Jansen      integer newflag, i, j
53059599516SKenneth E. Jansen      integer nsidl             ! number of surface IDs on-proc
53159599516SKenneth E. Jansen      integer nsidt             ! sum of all nsidl's
53259599516SKenneth E. Jansen      integer nsidg             ! number of surface IDs globally
53359599516SKenneth E. Jansen      integer nsid(numpe)       ! array of nsidl's
53459599516SKenneth E. Jansen      integer idisp(numpe)      ! needed by mpi: see note below
53559599516SKenneth E. Jansen      type llnod                ! structure for one node in a linked list
53659599516SKenneth E. Jansen        integer :: value
53759599516SKenneth E. Jansen        type (llnod), pointer :: next
53859599516SKenneth E. Jansen      end type llnod
53959599516SKenneth E. Jansen      type (llnod), pointer :: sidlist ! points to first elt of linked list
54059599516SKenneth E. Jansen      type (llnod), pointer :: sidelt  ! points to generic elt of linked list
5416f04c980SCameron Smith      type (llnod), pointer :: nextelt ! points to generic elt of linked list
54259599516SKenneth E. Jansen      integer, allocatable :: sidmapl(:) ! list of surfID's on-proc
54359599516SKenneth E. Jansen      integer, allocatable :: temp(:)    ! temp space
54459599516SKenneth E. Jansenc      integer iBCB(numelb,ndiBCB) ! column 1: naturalBC switches
54559599516SKenneth E. Jansen                                   ! column 2: surface ID's
54659599516SKenneth E. Jansenc Since we don't know a priori how many surface ID's there are,
54759599516SKenneth E. Jansenc on-processor or globally, we will store the ID's as an open-ended
54859599516SKenneth E. Jansenc link list while we determine the total number of distinct ID's
54959599516SKenneth E. Jansen      allocate (sidlist) ! allocate the first element of the sid
55059599516SKenneth E. Jansen                         ! linked list and point sidlist to it
55159599516SKenneth E. Jansen      nsidl=0            ! initialize # of sids on this processor
55259599516SKenneth E. Jansen      nsidg=0
5536f04c980SCameron Smith      nullify(sidlist % next)    ! next does not exist yet
55459599516SKenneth E. Jansen      do iblk=1, nelblb  ! loop over boundary element blocks
55559599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
55659599516SKenneth E. Jansen         do i = 1, npro         ! loop over boundary elements (belts)
55759599516SKenneth E. Jansen            iBCB2=miBCB(iblk)%p(i,2)
55859599516SKenneth E. Jansen            if(iBCB2.ne.zero) then ! if a sid is given for this belt
55959599516SKenneth E. Jansen               if(nsidl.eq.0) then     !   if this is the first sid we've seen
56059599516SKenneth E. Jansen                  nsidl=1              !       increment our count and
56159599516SKenneth E. Jansen                  sidlist % value = iBCB2    ! record its value
5626f04c980SCameron Smith                  nullify(sidlist % next)    ! next does not exist yet
56359599516SKenneth E. Jansen               else                    !   if this isn't the first sid
56459599516SKenneth E. Jansen                  newflag = 1          !     assume this is a new sid
56559599516SKenneth E. Jansen                  sidelt => sidlist    !     starting with the first sid
56659599516SKenneth E. Jansen                  do j = 1, nsidl      !     check the assumption
56759599516SKenneth E. Jansen                     if((sidelt % value).eq.iBCB2) then
56859599516SKenneth E. Jansen                        newflag=0      !        ...
56959599516SKenneth E. Jansen                     endif
57059599516SKenneth E. Jansen                     if(j.ne.nsidl) then !      ...
57159599516SKenneth E. Jansen                        sidelt => sidelt % next
57259599516SKenneth E. Jansen                     endif!                     ...
57359599516SKenneth E. Jansen                  enddo
57459599516SKenneth E. Jansen                  if(newflag.eq.1) then!     if it really is new to us
57559599516SKenneth E. Jansen                     nsidl = nsidl + 1 !         increment our count
57659599516SKenneth E. Jansen                     allocate (sidelt % next)!   tack a new element to the end
57759599516SKenneth E. Jansen                     sidelt => sidelt % next!    point to the new element
57859599516SKenneth E. Jansen                     sidelt % value = iBCB2    ! record the new sid
5796f04c980SCameron Smith                     nullify(sidelt % next)    ! next does not exist yet
58059599516SKenneth E. Jansen                  endif ! (really is a new sid)
58159599516SKenneth E. Jansen               endif ! (first sid)
58259599516SKenneth E. Jansen            endif ! (belt has a sid)
58359599516SKenneth E. Jansen         enddo ! (loop over belts)
58459599516SKenneth E. Jansen      enddo ! (loop over boundary element blocks)
58559599516SKenneth E. Jansenc Copy the data from the linked list to a more easily-used array form
58659599516SKenneth E. Jansen      if(nsidl.gt.0) then
58759599516SKenneth E. Jansen         allocate( sidmapl(nsidl) )
58859599516SKenneth E. Jansen         sidelt => sidlist      !     starting with the first sid
58959599516SKenneth E. Jansen         do j = 1, nsidl
59059599516SKenneth E. Jansen            sidmapl(j)=sidelt%value
59159599516SKenneth E. Jansen            if(j.ne.nsidl) sidelt => sidelt%next
59259599516SKenneth E. Jansen         enddo
59359599516SKenneth E. Jansen      else
59459599516SKenneth E. Jansen         allocate( sidmapl(1)) ! some compilers/MPI don't like to send unallocated arrays
59559599516SKenneth E. Jansen      endif
5966f04c980SCameron Smith!     Deallocate the link list
5976f04c980SCameron Smith!     http://stackoverflow.com/questions/9184675/how-does-fortran-deallocate-linked-lists
5986f04c980SCameron Smith      sidelt => sidlist
5996f04c980SCameron Smith      nextelt => sidelt % next
6006f04c980SCameron Smith      do
6016f04c980SCameron Smith        deallocate(sidelt)
6026f04c980SCameron Smith        if( .not. associated(nextelt) ) exit
6036f04c980SCameron Smith        sidelt => nextelt
6046f04c980SCameron Smith        nextelt => sidelt % next
6056f04c980SCameron Smith      enddo
6066f04c980SCameron Smith
60759599516SKenneth E. Jansenc Gather the number of surface ID's on each processor
60859599516SKenneth E. Jansen      if (numpe.gt.1) then      ! multiple processors
60959599516SKenneth E. Jansenc write the number of surfID's on the jth processor to slot j of nsid
61059599516SKenneth E. Jansen         call MPI_ALLGATHER(nsidl,1,MPI_INTEGER,nsid,1,
61159599516SKenneth E. Jansen     &          MPI_INTEGER,MPI_COMM_WORLD,ierr)
61259599516SKenneth E. Jansenc count up the total number of surfID's among all processes
61359599516SKenneth E. Jansen         nsidt=0
61459599516SKenneth E. Jansen         do j=1,numpe
61559599516SKenneth E. Jansen            nsidt=nsidt+nsid(j)
61659599516SKenneth E. Jansen         enddo
61759599516SKenneth E. Jansen      else                      ! single processor
61859599516SKenneth E. Jansenc the local information is the global information for single-processor
61959599516SKenneth E. Jansen         nsid=nsidl
62059599516SKenneth E. Jansen         nsidt=nsidl
62159599516SKenneth E. Jansen      endif                     ! if-else for multiple processors
62259599516SKenneth E. Jansen      if(nsidt.gt.0) then
62359599516SKenneth E. Jansenc
62459599516SKenneth E. Jansenc  Make all-processor surfID collage
62559599516SKenneth E. Jansenc
62659599516SKenneth E. Jansenc there will be some duplicate surface ID's when we gather, so we
62759599516SKenneth E. Jansenc will use a temporary array
62859599516SKenneth E. Jansen         allocate (temp(nsidt))
62959599516SKenneth E. Jansen         if (numpe.gt.1) then   ! multiple processors
63059599516SKenneth E. Jansenc we will gather surfID's from local on-proc sets to a global set
63159599516SKenneth E. Jansenc we will stack each processor's surfID list atop that of the previous
63259599516SKenneth E. Jansenc processor.  If the list for processor i is called sidmapi, then our
63359599516SKenneth E. Jansenc global coordinate list sidmap will look like this:
63459599516SKenneth E. Jansenc ---------------------------------------------------------------
63559599516SKenneth E. Jansenc | sidmap1       | sidmap2           | sidmap3   |   ...       |
63659599516SKenneth E. Jansenc ---------------------------------------------------------------
63759599516SKenneth E. Jansenc  <---nsid(1)---> <-----nsid(2)-----> <-nsid(3)->
63859599516SKenneth E. Jansenc  <------------------------nsidt-----------------------...---->
63959599516SKenneth E. Jansenc To accomplish this with MPI, we use MPI_ALLGATHERV, summarized as:
64059599516SKenneth E. JansencMPI_ALLGATHERV(sendbuf,sendcount,sendtype,recvbuf,recvcount,disp,recvtype,comm)
64159599516SKenneth E. Jansenc[ IN sendbuf] starting address of send buffer (choice)
64259599516SKenneth E. Jansenc[ IN sendcount] number of elements in send buffer (integer)
64359599516SKenneth E. Jansenc[ IN sendtype] data type of send buffer elements (handle)
64459599516SKenneth E. Jansenc[ OUT recvbuf] address of receive buffer (choice)
64559599516SKenneth E. Jansenc[ IN recvcount] number of elements received from each process (int array)
64659599516SKenneth E. Jansenc[ IN disp] displacement array
64759599516SKenneth E. Jansenc[ IN recvtype] data type of receive buffer elements (handle)
64859599516SKenneth E. Jansenc[ IN comm] communicator (handle)
64959599516SKenneth E. Jansenc The displacement array disp is an array of integers in which the jth
65059599516SKenneth E. Jansenc entry indicates which slot of sidmap marks beginning of sidmapj
65159599516SKenneth E. Jansenc So, first we will build this displacement array
65259599516SKenneth E. Jansen            idisp(:)=0      ! starting with zero, since MPI likes C-numbering
65359599516SKenneth E. Jansen            do j=2,numpe
65459599516SKenneth E. Jansen               idisp(j)=idisp(j-1)+nsid(j-1) ! see diagram above
65559599516SKenneth E. Jansen            enddo
65659599516SKenneth E. Jansenc Now, we gather the data
65759599516SKenneth E. Jansen            call MPI_ALLGATHERV(sidmapl(:),nsidl,
65859599516SKenneth E. Jansen     &           MPI_INTEGER,temp(:),nsid,idisp,
65959599516SKenneth E. Jansen     &           MPI_INTEGER,MPI_COMM_WORLD,ierr)
66059599516SKenneth E. Jansenc sort surfID's, lowest to highest
66159599516SKenneth E. Jansen            isorted = 0
66259599516SKenneth E. Jansen            do while (isorted.eq.0) ! while the list isn't sorted
66359599516SKenneth E. Jansen               isorted = 1      ! assume the list is sorted this time
66459599516SKenneth E. Jansen               do j = 2, nsidt  ! loop over ID's
66559599516SKenneth E. Jansen                  if(temp(j).lt.temp(j-1)) then ! ID exceeds predecessor
66659599516SKenneth E. Jansen                     itmp=temp(j-1)
66759599516SKenneth E. Jansen                     temp(j-1)=temp(j)
66859599516SKenneth E. Jansen                     temp(j)=itmp
66959599516SKenneth E. Jansen                     isorted=0  ! not sorted yet
67059599516SKenneth E. Jansen                  endif
67159599516SKenneth E. Jansen               enddo            !loop over ID's
67259599516SKenneth E. Jansen            enddo               ! while not sorted
67359599516SKenneth E. Jansenc determine the total number of distinct surfID's, globally
67459599516SKenneth E. Jansen            nsidg=nsidt         ! assume there are no duplicate SurfID's
67559599516SKenneth E. Jansen            do j = 2, nsidt
67659599516SKenneth E. Jansen               if(temp(j).eq.temp(j-1)) nsidg = nsidg - 1 ! correction
67759599516SKenneth E. Jansen            enddo
67859599516SKenneth E. Jansenc create duplicate-free surfID list
67959599516SKenneth E. Jansen            allocate( sidmapg(nsidg) )
68059599516SKenneth E. Jansen            sidmapg(1)=temp(1)
68159599516SKenneth E. Jansen            nsidg = 1
68259599516SKenneth E. Jansen            do j = 2, nsidt
68359599516SKenneth E. Jansen               if(temp(j).ne.temp(j-1)) then
68459599516SKenneth E. Jansen                  nsidg = nsidg + 1
68559599516SKenneth E. Jansen                  sidmapg(nsidg)=temp(j)
68659599516SKenneth E. Jansen               endif
68759599516SKenneth E. Jansen            enddo
68859599516SKenneth E. Jansen            deallocate( temp )
68959599516SKenneth E. Jansen         else                   ! single-processor
69059599516SKenneth E. Jansenc global data is local data in single processor case
69159599516SKenneth E. Jansen            nsidg=nsidl
69259599516SKenneth E. Jansen            allocate( sidmapg(nsidg) )
69359599516SKenneth E. Jansen            sidmapg=sidmapl
6946f04c980SCameron Smith            deallocate(sidmapl)
69559599516SKenneth E. Jansen         endif                  ! if-else multiple processor
69659599516SKenneth E. Jansenc
69759599516SKenneth E. Jansen      endif                     ! if at least one surfid
69859599516SKenneth E. Jansenc
69959599516SKenneth E. Jansen      return
70059599516SKenneth E. Jansenc
70159599516SKenneth E. Jansen      end
70259599516SKenneth E. Jansen
70359599516SKenneth E. Jansen
70459599516SKenneth E. Jansen
70559599516SKenneth E. Jansen      subroutine genotwn(x, BCtmp, iBC, nsurf)
70659599516SKenneth E. Jansenc
70759599516SKenneth E. Jansenc----------------------------------------------------------------------
70859599516SKenneth E. Jansenc  This routine determines which node to use as the first node off the
70959599516SKenneth E. Jansenc  wall in near-wall modeling traction calculations for each wall node.
71059599516SKenneth E. Jansenc  Each wall node has a normal, calculated from the wall elements to
71159599516SKenneth E. Jansenc  which that node belongs.  We seek the node within the boundary
71259599516SKenneth E. Jansenc  element that lies closest to the line defined by the normal vector.
71359599516SKenneth E. Jansenc  We create normalized vectors pointing from the wall node in
71459599516SKenneth E. Jansenc  question to each of the nodes in the boundary element. The vector
71559599516SKenneth E. Jansenc  that has the largest projection onto the wall node's normal vector
71659599516SKenneth E. Jansenc  points to the node we want.  Nodes that are not on a wall point to
71759599516SKenneth E. Jansenc  themselves as their own off-the-wall node.
71859599516SKenneth E. Jansenc
71959599516SKenneth E. Jansenc input:
72059599516SKenneth E. Jansenc  x     (nshg,3)     :        : nodal position vectors
72159599516SKenneth E. Jansenc  wnrm  (nshg,3)     : (mod)  : normal vectors for each node
72259599516SKenneth E. Jansenc  iBCB5 (nshg)       : (file) : wall flag for belt
72359599516SKenneth E. Jansenc  ienb  (numelb,nen) : (file) : connectivity for belts
72459599516SKenneth E. Jansenc
72559599516SKenneth E. Jansenc output:
72659599516SKenneth E. Jansenc  otwn  (nshg)       : (mod)  : off-the-wall nodes for each node
72759599516SKenneth E. Jansenc
72859599516SKenneth E. Jansenc
72959599516SKenneth E. Jansenc Kenneth Jansen, Summer 2000 (algorithm)
73059599516SKenneth E. Jansenc Michael Yaworski, Summer 2000 (code)
73159599516SKenneth E. Jansenc----------------------------------------------------------------------
73259599516SKenneth E. Jansenc
73359599516SKenneth E. Jansen      use pointer_data          ! used for mienb, miBCB
73459599516SKenneth E. Jansen      use turbSA
73559599516SKenneth E. Jansen      include "common.h"
73659599516SKenneth E. Jansenc
73759599516SKenneth E. Jansen      integer iel, nod, can
73859599516SKenneth E. Jansen      real*8 vec(3), leng, dp, bigdp, lil
73959599516SKenneth E. Jansen      real*8 x(numnp,nsd),BCtmp(nshg,ndof+7)
74059599516SKenneth E. Jansen      integer iBC(nshg), nsurf(nshg)
74159599516SKenneth E. Jansen      integer gbits
74259599516SKenneth E. Jansen      integer, allocatable :: ienb(:)
74359599516SKenneth E. Jansen
74459599516SKenneth E. Jansen      allocate( otwn(nshg) )
74559599516SKenneth E. Jansenc
74659599516SKenneth E. Jansenc initialize otwn to oneself
74759599516SKenneth E. Jansenc
74859599516SKenneth E. Jansen      do nod = 1, nshg
74959599516SKenneth E. Jansen         otwn(nod)=nod
75059599516SKenneth E. Jansen      enddo
75159599516SKenneth E. Jansenc
75259599516SKenneth E. Jansenc determine otwn for each wall node
75359599516SKenneth E. Jansenc
75459599516SKenneth E. Jansen      do iblk=1, nelblb         ! loop over boundary element blocks
75559599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
75659599516SKenneth E. Jansen         nenl  = lcblkb(5,iblk)
75759599516SKenneth E. Jansen         nenbl = lcblkb(6,iblk)
75859599516SKenneth E. Jansen         nshl = lcblkb(9,iblk)
75959599516SKenneth E. Jansen         allocate( ienb(nshl) )
76059599516SKenneth E. Jansen         do i = 1, npro         ! loop over boundary elements
76159599516SKenneth E. Jansen            iBCB1=miBCB(iblk)%p(i,1)
76259599516SKenneth E. Jansen            iBCB2=miBCB(iblk)%p(i,2)
76359599516SKenneth E. Jansen            ienb(1:nshl)=mienb(iblk)%p(i,1:nshl)
76459599516SKenneth E. Jansen            if (btest(iBCB1,4)) then ! (on the wall)
76559599516SKenneth E. Jansen               do nod = 1, nenbl !   for each wall node in this belt
76659599516SKenneth E. Jansen                  bigdp = zero  !     initialize largest dot product
76759599516SKenneth E. Jansen                  do can=nenbl+1,nenl !  loop over off-the-wall candidates
76859599516SKenneth E. Jansen                     nn=ienb(can)
76959599516SKenneth E. Jansen                               !       don't bother with wall nodes
77059599516SKenneth E. Jansen                     if(nsurf(nn).ne.0) cycle
77159599516SKenneth E. Jansen                               !       don't bother with no-slip nodes
77259599516SKenneth E. Jansen                     if(ibits(iBC(nn),3,3).eq.7 .and.
77359599516SKenneth E. Jansen     &                    BCtmp(nn,7).eq.zero) cycle
77459599516SKenneth E. Jansenc                              !       calculate candidate vector
77559599516SKenneth E. Jansen                     vec(:)=x(ienb(can),:)-x(ienb(nod),:)
77659599516SKenneth E. Jansen                     leng=sqrt(vec(1)**2+vec(2)**2+vec(3)**2)
77759599516SKenneth E. Jansen                     vec(:)=vec(:)/leng
77859599516SKenneth E. Jansenc                              !       calculate dot product with wnrm
77959599516SKenneth E. Jansen                     vec(:)=vec(:)*wnrm(ienb(nod),:)
78059599516SKenneth E. Jansen                     dp=vec(1)+vec(2)+vec(3)
78159599516SKenneth E. Jansenc                              !       set candidate as otwn if necessary
78259599516SKenneth E. Jansenc                              !       (wnrm points into fluid, so
78359599516SKenneth E. Jansenc                              !        we want the most positive dp)
78459599516SKenneth E. Jansen                     if (dp.gt.bigdp) then
78559599516SKenneth E. Jansen                        otwn(ienb(nod)) = ienb(can)
78659599516SKenneth E. Jansen                        bigdp=dp
78759599516SKenneth E. Jansen                     endif
78859599516SKenneth E. Jansen                  enddo         !(loop over off-the-wall candidates)
78959599516SKenneth E. Jansen               enddo            !(loop over wall nodes in current belt)
79059599516SKenneth E. Jansen            endif
79159599516SKenneth E. Jansen         enddo                  !(loop over elts in block)
79259599516SKenneth E. Jansen         deallocate(ienb)
79359599516SKenneth E. Jansen      enddo                     !(loop over boundary element blocks)
79459599516SKenneth E. Jansen      do nn = 1, nshg
79559599516SKenneth E. Jansen         if((otwn(nn).eq.nn).and.(nsurf(nn).ne.0)) then ! if a node on a
79659599516SKenneth E. Jansen                                                    ! modeled surface
79759599516SKenneth E. Jansen                                                    ! didn't find a node
79859599516SKenneth E. Jansen                                                    ! off the wall
79959599516SKenneth E. Jansen            lil = 1.0e32 ! distance to current closest prospect
80059599516SKenneth E. Jansen            do can = 1, nshg ! loop over nodes
80159599516SKenneth E. Jansen               if(nsurf(can).eq.0) then ! if this candidate is off the
80259599516SKenneth E. Jansen                                        ! wall
80359599516SKenneth E. Jansen                  if(ibits(iBC(can),3,3).eq.7 .and.
80459599516SKenneth E. Jansen     &               BCtmp(can,7).eq.zero) then  ! no slip nodes not allowed
80559599516SKenneth E. Jansen                  else          ! not a forbidden node
80659599516SKenneth E. Jansen                     vec(:)=x(nn,:)-x(can,:)
80759599516SKenneth E. Jansen                     leng=vec(1)**2+vec(2)**2+vec(3)**2
80859599516SKenneth E. Jansen                     if(leng.lt.lil) then ! this candidate is closest so far
80959599516SKenneth E. Jansen                        lil=leng
81059599516SKenneth E. Jansen                        otwn(nn)=can
81159599516SKenneth E. Jansen                     endif      ! closest so far
81259599516SKenneth E. Jansen                  endif  ! end of no slip nodes
81359599516SKenneth E. Jansen               endif ! off-the-wall check
81459599516SKenneth E. Jansen            enddo ! loop over nodes
81559599516SKenneth E. Jansen         endif ! lonely wall-node check
81659599516SKenneth E. Jansen      enddo ! loop over nodes
81759599516SKenneth E. Jansenc
81859599516SKenneth E. Jansen      return
81959599516SKenneth E. Jansenc
82059599516SKenneth E. Jansen      end
82159599516SKenneth E. Jansen
822