xref: /phasta/phSolver/common/genbc.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen      subroutine genBC (iBC,  BC,   x,   ilwork, iper)
2*59599516SKenneth E. Jansenc
3*59599516SKenneth E. Jansenc----------------------------------------------------------------------
4*59599516SKenneth E. Jansenc  This routine generates the essential prescribed boundary conditions.
5*59599516SKenneth E. Jansenc
6*59599516SKenneth E. Jansenc input:
7*59599516SKenneth E. Jansenc  iBC   (nshg)        : boundary condition code
8*59599516SKenneth E. Jansenc  nBC   (nshg)        : boundary condition mapping array
9*59599516SKenneth E. Jansenc
10*59599516SKenneth E. Jansenc output:
11*59599516SKenneth E. Jansenc  BC    (nshg,ndofBC) : The constraint data for prescribed BC
12*59599516SKenneth E. Jansenc
13*59599516SKenneth E. Jansenc
14*59599516SKenneth E. Jansenc Note: genBC1 reduces the input data for the velocity. In the
15*59599516SKenneth E. Jansenc       case of varying velocity direction in the generation, the
16*59599516SKenneth E. Jansenc       results may not be correct. (since a linearity assumption is
17*59599516SKenneth E. Jansenc       made in the generation).
18*59599516SKenneth E. Jansenc
19*59599516SKenneth E. Jansenc
20*59599516SKenneth E. Jansenc Farzin Shakib, Spring 1986.
21*59599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
22*59599516SKenneth E. Jansenc----------------------------------------------------------------------
23*59599516SKenneth E. Jansenc
24*59599516SKenneth E. Jansen      use slpw
25*59599516SKenneth E. Jansen      use readarrays            ! used to access BCinp, nBC
26*59599516SKenneth E. Jansen      use specialBC ! filling acs here
27*59599516SKenneth E. Jansen      include "common.h"
28*59599516SKenneth E. Jansenc
29*59599516SKenneth E. Jansen      dimension iBC(nshg),                nsurf(nshg),
30*59599516SKenneth E. Jansen     &            BC(nshg,ndofBC),
31*59599516SKenneth E. Jansen     &            x(numnp,nsd),             ilwork(nlwork),
32*59599516SKenneth E. Jansen     &            iper(nshg)
33*59599516SKenneth E. Jansenc
34*59599516SKenneth E. Jansenc BCinp for each point has:
35*59599516SKenneth E. Jansenc   D T P c11 c12 c13 M1 c21 c22 c23 M2 theta S1 S2 S3...
36*59599516SKenneth E. Jansenc   1 2 3 4   5   6   7  8   9   10  11 12    13 14 15...
37*59599516SKenneth E. Jansenc Remember, ndof=nsd+2+nsclr
38*59599516SKenneth E. Jansenc
39*59599516SKenneth E. Jansenc Arrays in the following 1 line are now dimensioned in readnblk
40*59599516SKenneth E. Jansenc        dimension BCinp(numpbc,ndof+7)
41*59599516SKenneth E. Jansenc
42*59599516SKenneth E. Jansen      dimension BCtmp(nshg,ndof+7)
43*59599516SKenneth E. Jansenc
44*59599516SKenneth E. Jansenc ndof+7= 3(thermos) + (nsd-1)*(nsd+1) + nscalars + 1 (theta)
45*59599516SKenneth E. Jansenc                       #vect *(vec dir +mag)
46*59599516SKenneth E. Jansenc
47*59599516SKenneth E. Jansenc.... --------------------------->  Input  <---------------------------
48*59599516SKenneth E. Jansenc
49*59599516SKenneth E. Jansenc.... convert boundary condition data
50*59599516SKenneth E. Jansenc
51*59599516SKenneth E. Jansen      BCtmp = zero
52*59599516SKenneth E. Jansenc
53*59599516SKenneth E. Jansen      if(numpbc.ne.0) then
54*59599516SKenneth E. Jansen         do i = 1, ndof+7
55*59599516SKenneth E. Jansen            where (nBC(:) .ne. 0) BCtmp(:,i) = BCinp(nBC(:),i)
56*59599516SKenneth E. Jansen         enddo
57*59599516SKenneth E. Jansen         deallocate(BCinp)
58*59599516SKenneth E. Jansen      endif
59*59599516SKenneth E. Jansen
60*59599516SKenneth E. Jansenc
61*59599516SKenneth E. Jansen      if(any(BCtmp(:,12).ne.0)) then
62*59599516SKenneth E. Jansen         iabc=1
63*59599516SKenneth E. Jansen         allocate (acs(nshg,2))
64*59599516SKenneth E. Jansen         where (btest(iBC,10))
65*59599516SKenneth E. Jansen            acs(:,1) = cos(BCtmp(:,12))
66*59599516SKenneth E. Jansen            acs(:,2) = sin(BCtmp(:,12))
67*59599516SKenneth E. Jansen         endwhere
68*59599516SKenneth E. Jansen      endif
69*59599516SKenneth E. Jansen
70*59599516SKenneth E. Jansenc
71*59599516SKenneth E. Jansenc.... ----------------------> Wall Normals  <--------------------------
72*59599516SKenneth E. Jansenc (calculate the normal and adjust BCinp to the true normal as needed)
73*59599516SKenneth E. Jansenc
74*59599516SKenneth E. Jansen      if(navier.eq.1)then
75*59599516SKenneth E. Jansen         call genwnm (iBC,  BCtmp,   x,   ilwork, iper, nsurf)
76*59599516SKenneth E. Jansen      endif
77*59599516SKenneth E. Jansenc  determine the first point off the wall for each wall node
78*59599516SKenneth E. Jansen      if(navier.eq.1)then
79*59599516SKenneth E. Jansen         call genotwn (x,BCtmp, iBC, nsurf)
80*59599516SKenneth E. Jansen      endif
81*59599516SKenneth E. Jansenc.... ------------------------>  Conversion  <-------------------------
82*59599516SKenneth E. Jansenc
83*59599516SKenneth E. Jansenc.... convert the input boundary conditions to condensed version
84*59599516SKenneth E. Jansenc
85*59599516SKenneth E. Jansen      BC = zero
86*59599516SKenneth E. Jansenc
87*59599516SKenneth E. Jansen      if(navier.eq.1)then ! zero navier means Euler simulation
88*59599516SKenneth E. Jansen         call genBC1 (BCtmp,  iBC,  BC)
89*59599516SKenneth E. Jansen      else ! navier equals zero
90*59599516SKenneth E. Jansen         allocate(BCtmpg(nshg,ndof+7))
91*59599516SKenneth E. Jansen         allocate(BCg(nshg,ndofBC))
92*59599516SKenneth E. Jansen         allocate(iBCg(nshg))
93*59599516SKenneth E. Jansen         BCtmpg=BCtmp
94*59599516SKenneth E. Jansen         iBCg=iBC
95*59599516SKenneth E. Jansen         call genBC1 (BCtmp,  iBC,  BC)
96*59599516SKenneth E. Jansenc... genBC1 convert BCtmp to BC
97*59599516SKenneth E. Jansen         BCg=BC
98*59599516SKenneth E. Jansen         icdg=icd
99*59599516SKenneth E. Jansenc... find slip wall
100*59599516SKenneth E. Jansen         call findslpw(x,ilwork,iper,iBC)
101*59599516SKenneth E. Jansenc... apply slip wall condition to wall nodes
102*59599516SKenneth E. Jansen         do i=1,nslpwnd
103*59599516SKenneth E. Jansen            nn=idxslpw(i)
104*59599516SKenneth E. Jansen            call slpwBC(mpslpw(nn,1),mpslpw(nn,2),iBCg(nn),
105*59599516SKenneth E. Jansen     &               BCg(nn,:),  BCtmpg(nn,:),
106*59599516SKenneth E. Jansen     &               iBC(nn),    BC(nn,:),
107*59599516SKenneth E. Jansen     &               wlnorm(nn,:,:)                     )
108*59599516SKenneth E. Jansen         enddo
109*59599516SKenneth E. Jansen         icd=icdg
110*59599516SKenneth E. Jansen      endif
111*59599516SKenneth E. Jansen
112*59599516SKenneth E. Jansen
113*59599516SKenneth E. Jansenc
114*59599516SKenneth E. Jansenc.... --------------------------->  Echo  <----------------------------
115*59599516SKenneth E. Jansenc
116*59599516SKenneth E. Jansenc.... echo the input data
117*59599516SKenneth E. Jansenc
118*59599516SKenneth E. Jansen      if (necho .lt. 3) then
119*59599516SKenneth E. Jansen         nn = 0
120*59599516SKenneth E. Jansen         do n = 1, nshg
121*59599516SKenneth E. Jansen            if (nBC(n) .ne. 0) then
122*59599516SKenneth E. Jansen               nn = nn + 1
123*59599516SKenneth E. Jansen               if(mod(nn,50).eq.1)
124*59599516SKenneth E. Jansen     &              write(iecho,1000)ititle,(j,j=1,ndofBC)
125*59599516SKenneth E. Jansen               write (iecho,1100) n, (BC(n,i),i=1,ndofBC)
126*59599516SKenneth E. Jansen            endif
127*59599516SKenneth E. Jansen         enddo
128*59599516SKenneth E. Jansen      endif
129*59599516SKenneth E. Jansenc
130*59599516SKenneth E. Jansenc.... return
131*59599516SKenneth E. Jansenc
132*59599516SKenneth E. Jansen      return
133*59599516SKenneth E. Jansenc
134*59599516SKenneth E. Jansen 1000 format(a80,//,
135*59599516SKenneth 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',//,
136*59599516SKenneth E. Jansen     &  '    Node  ',/,
137*59599516SKenneth E. Jansen     &  '   Number ',5x,6('BC',i1,:,10x))
138*59599516SKenneth E. Jansen 1100 format(1p,2x,i5,3x,6(e12.5,1x))
139*59599516SKenneth E. Jansenc
140*59599516SKenneth E. Jansen      end
141*59599516SKenneth E. Jansen
142*59599516SKenneth E. Jansen
143*59599516SKenneth E. Jansen
144*59599516SKenneth E. Jansen
145*59599516SKenneth E. Jansen
146*59599516SKenneth E. Jansen
147*59599516SKenneth E. Jansen
148*59599516SKenneth E. Jansen
149*59599516SKenneth E. Jansen
150*59599516SKenneth E. Jansen
151*59599516SKenneth E. Jansen      subroutine genwnm (iBC, BCtmp,    x,   ilwork, iper, nsurf)
152*59599516SKenneth E. Jansenc----------------------------------------------------------------------
153*59599516SKenneth E. Jansenc  This routine generates the normal to a wall
154*59599516SKenneth E. Jansenc
155*59599516SKenneth E. Jansenc input:
156*59599516SKenneth E. Jansenc  iBC   (nshg)        : boundary condition code
157*59599516SKenneth E. Jansenc  nBC   (nshg)        : boundary condition mapping array
158*59599516SKenneth E. Jansenc
159*59599516SKenneth E. Jansenc output:
160*59599516SKenneth E. Jansenc  BCtmp    (nshg,ndof+6) : The constraint data for prescribed BC
161*59599516SKenneth E. Jansenc
162*59599516SKenneth E. Jansenc
163*59599516SKenneth E. Jansenc----------------------------------------------------------------------
164*59599516SKenneth E. Jansenc
165*59599516SKenneth E. Jansen      use turbSA
166*59599516SKenneth E. Jansen      use pointer_data          ! used for mienb, mibcb
167*59599516SKenneth E. Jansen      include "common.h"
168*59599516SKenneth E. Jansen      include "mpif.h"
169*59599516SKenneth E. Jansenc
170*59599516SKenneth E. Jansen      character*20 fname1,  fmt1
171*59599516SKenneth E. Jansen      character*5  cname
172*59599516SKenneth E. Jansen      dimension iBC(nshg),                  iper(nshg),
173*59599516SKenneth E. Jansen     &            x(numnp,nsd),             ilwork(nlwork)
174*59599516SKenneth E. Jansenc
175*59599516SKenneth E. Jansen      dimension  BCtmpSAV(nshg,ndof+7)
176*59599516SKenneth E. Jansen      dimension  BCtmp(nshg,ndof+7),      fBC(nshg,ndofBC),
177*59599516SKenneth E. Jansen     &            e1(3),                    e2(3),
178*59599516SKenneth E. Jansen     &            elnrm(3),                 asum(numnp)
179*59599516SKenneth E. Jansenc
180*59599516SKenneth E. Jansen      integer sid, nsidg
181*59599516SKenneth E. Jansen      integer nsurf(nshg), ivec(nsd)
182*59599516SKenneth E. Jansen      logical :: firstvisit(nshg)
183*59599516SKenneth E. Jansen      real*8 BCvecs(2,nsd)
184*59599516SKenneth E. Jansen      integer, allocatable :: ienb(:)
185*59599516SKenneth E. Jansen      dimension wnmdb(nshg,nsd)
186*59599516SKenneth E. Jansenc
187*59599516SKenneth E. Jansenc  wnrm is dimensioned nshg but the math is only done for straight sided
188*59599516SKenneth E. Jansenc  elements at this point so wnrm will not be calculated for the hierarchic
189*59599516SKenneth E. Jansenc  modes.  Note that the wall model creates a p.w. linear representation
190*59599516SKenneth E. Jansenc  only at this time.
191*59599516SKenneth E. Jansenc
192*59599516SKenneth E. Jansen      allocate ( wnrm(nshg,3) )
193*59599516SKenneth E. Jansenc
194*59599516SKenneth E. Jansenc.... ----------------------> Wall Normals  <--------------------------
195*59599516SKenneth E. Jansenc (calculate the normal and adjust BCinp to the true normal as needed)
196*59599516SKenneth E. Jansenc
197*59599516SKenneth E. Jansenc
198*59599516SKenneth E. Jansen      asum = zero
199*59599516SKenneth E. Jansen      wnrm = zero
200*59599516SKenneth E. Jansen
201*59599516SKenneth E. Jansenc
202*59599516SKenneth E. Jansenc....  Save a copy of BCtmp so that after we calculate the normals we
203*59599516SKenneth E. Jansenc      can recover the comp3 information.
204*59599516SKenneth E. Jansenc
205*59599516SKenneth E. Jansen      BCtmpSAV=BCtmp
206*59599516SKenneth E. Jansenc
207*59599516SKenneth E. Jansenc Count out the number of surface ID's on-processor, and map them
208*59599516SKenneth E. Jansen      call gensidcount(nsidg)
209*59599516SKenneth E. Jansen
210*59599516SKenneth E. Jansenc
211*59599516SKenneth E. Jansen      if(nsidg.gt.0) then       ! if there are any surfID's
212*59599516SKenneth E. Jansen         nsurf(:) = 0
213*59599516SKenneth E. Jansen         do k = 1, nsidg        ! loop over Surface ID's
214*59599516SKenneth E. Jansen            sid = sidmapg(k)
215*59599516SKenneth E. Jansen            firstvisit(:)=.true.
216*59599516SKenneth E. Jansen            wnrm(:,1:3)=zero
217*59599516SKenneth E. Jansen            do iblk=1, nelblb   ! loop over boundary element blocks
218*59599516SKenneth E. Jansen               npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
219*59599516SKenneth E. Jansen               nenbl = lcblkb(6,iblk)
220*59599516SKenneth E. Jansen               nshl = lcblkb(9,iblk)
221*59599516SKenneth E. Jansen               allocate( ienb(nshl) )
222*59599516SKenneth E. Jansen               do i = 1, npro   ! loop over boundary elements
223*59599516SKenneth E. Jansen                  iBCB1=miBCB(iblk)%p(i,1)
224*59599516SKenneth E. Jansen                  iBCB2=miBCB(iblk)%p(i,2)
225*59599516SKenneth E. Jansen                  ienb(1:nshl)=mienb(iblk)%p(i,1:nshl)
226*59599516SKenneth E. Jansenc don't bother with elements that aren't on modeled surfaces
227*59599516SKenneth E. Jansenc              if ( not          (wall set  and   traction set)    )
228*59599516SKenneth E. Jansen                  if (.not.(btest(iBCB1,2).and.btest(iBCB1,4)))
229*59599516SKenneth E. Jansen     &                 cycle
230*59599516SKenneth E. Jansenc don't bother with elements that don't lie on the current surface
231*59599516SKenneth E. Jansen                  if (iBCB2.ne.sid) cycle
232*59599516SKenneth E. Jansenc
233*59599516SKenneth E. Jansenc.... calculate this element's area-weighted normal vector
234*59599516SKenneth E. Jansenc
235*59599516SKenneth E. Jansen                  e1 = x(ienb(2),:)-x(ienb(1),:)
236*59599516SKenneth E. Jansen                  e2 = x(ienb(3),:)-x(ienb(1),:)
237*59599516SKenneth E. Jansen                  elnrm(1) = e1(2)*e2(3)-e1(3)*e2(2)
238*59599516SKenneth E. Jansen                  elnrm(2) = e1(3)*e2(1)-e1(1)*e2(3)
239*59599516SKenneth E. Jansen                  elnrm(3) = e1(1)*e2(2)-e1(2)*e2(1)
240*59599516SKenneth E. Jansenc Tetrahedral elements have negative volumes in phastaI, so
241*59599516SKenneth E. Jansenc the normals calculated from the boundary faces must be inverted
242*59599516SKenneth E. Jansenc to point into the fluid
243*59599516SKenneth E. Jansen                  if(nenbl.eq.3) elnrm(:)=-elnrm(:)
244*59599516SKenneth E. Jansenc
245*59599516SKenneth E. Jansenc.... add area-weighted normals to the nodal tallies for this surface
246*59599516SKenneth E. Jansenc
247*59599516SKenneth E. Jansen                  do j = 1, nenbl ! loop over elt boundary nodes
248*59599516SKenneth E. Jansen                     nn=ienb(j) ! global node number
249*59599516SKenneth E. Jansen                     if(firstvisit(nn)) then
250*59599516SKenneth E. Jansen                        firstvisit(nn)=.false.
251*59599516SKenneth E. Jansen                        nsurf(nn)=nsurf(nn)+1
252*59599516SKenneth E. Jansen                        if(nsurf(nn).eq.1) BCtmp(nn,4:6)=zero
253*59599516SKenneth E. Jansen                        if(nsurf(nn).eq.2) BCtmp(nn,8:10)=zero
254*59599516SKenneth E. Jansen                     endif
255*59599516SKenneth E. Jansen                     wnrm(nn,:)=wnrm(nn,:)+elnrm
256*59599516SKenneth E. Jansen                  enddo         ! loop over elt boundary nodes
257*59599516SKenneth E. Jansen               enddo            ! end loop over boundary elements in block
258*59599516SKenneth E. Jansen               deallocate(ienb)
259*59599516SKenneth E. Jansen            enddo               ! end loop over boundary element blocks
260*59599516SKenneth E. Jansenc Now we have all of this surface's contributions to wall normals
261*59599516SKenneth E. Jansenc for all nodes, along with an indication of how many surfaces
262*59599516SKenneth E. Jansenc each node has encountered so far.  Contributions from other processors
263*59599516SKenneth E. Jansenc should now be accumulated for this surface
264*59599516SKenneth E. Jansenc
265*59599516SKenneth E. Jansenc axisymm BC's need BC (and we have not built it yet) so we need to create
266*59599516SKenneth E. Jansenc the entries it needs.
267*59599516SKenneth E. Jansenc
268*59599516SKenneth E. Jansen            if(iabc==1) then
269*59599516SKenneth E. Jansen               where (btest(iBC,10))
270*59599516SKenneth E. Jansen                  fBC(:,1) = cos(BCtmp(:,12))
271*59599516SKenneth E. Jansen                  fBC(:,2) = sin(BCtmp(:,12))
272*59599516SKenneth E. Jansen               endwhere
273*59599516SKenneth E. Jansen
274*59599516SKenneth E. Jansenc before the commu we need to rotate the residual vector for axisymmetric
275*59599516SKenneth E. Jansenc boundary conditions (so that off processor periodicity is a dof add instead
276*59599516SKenneth E. Jansenc of a dof combination).  Take care of all nodes now so periodicity, like
277*59599516SKenneth E. Jansenc commu is a simple dof add.
278*59599516SKenneth E. Jansen               call rotabc(wnrm, iBC, 'in ')
279*59599516SKenneth E. Jansen            endif
280*59599516SKenneth E. Jansen
281*59599516SKenneth E. Jansen            if (numpe.gt.1) then
282*59599516SKenneth E. Jansenc.... add areas and norms contributed from other processors
283*59599516SKenneth E. Jansen               call commu (wnrm, ilwork, 3, 'in ')
284*59599516SKenneth E. Jansen            endif
285*59599516SKenneth E. Jansenc.... account for periodicity and axisymmetry
286*59599516SKenneth E. Jansen            call bc3per(iBC,wnrm, iper,ilwork, 3)
287*59599516SKenneth E. Jansenc at this point the master has all the information (slaves are zeroed and
288*59599516SKenneth E. Jansenc waiting to be told what the master has...lets do that).
289*59599516SKenneth E. Jansenc
290*59599516SKenneth E. Jansenc.... local periodic (and axisymmetric) boundary conditions (no communications)
291*59599516SKenneth E. Jansen            wnmdb=wnrm
292*59599516SKenneth E. Jansen            do i = 1,numnp      ! only use the vertices in the normal calc
293*59599516SKenneth E. Jansen               wnrm(i,:) = wnrm(iper(i),:)
294*59599516SKenneth E. Jansen            enddo
295*59599516SKenneth E. Jansen            wnmdb=wnrm
296*59599516SKenneth E. Jansen            if (numpe.gt.1) then
297*59599516SKenneth E. Jansenc.... tell other processors the resulting (and now complete) sums
298*59599516SKenneth E. Jansen               call commu (wnrm, ilwork, 3, 'out')
299*59599516SKenneth E. Jansen            endif
300*59599516SKenneth E. Jansenc Axisymmetric slaves need to be rotated back
301*59599516SKenneth E. Jansen            if(iabc==1) then    !are there any axisym bc's
302*59599516SKenneth E. Jansen               call rotabc(wnrm, iBC, 'out')
303*59599516SKenneth E. Jansen            endif
304*59599516SKenneth E. Jansenc Now all nodes have all the normal contributions for this surface,
305*59599516SKenneth E. Jansenc including those from off-processor and periodic nodes.
306*59599516SKenneth E. Jansenc We distribute these summed vectors to the proper slots in BCtmp,
307*59599516SKenneth E. Jansenc according to how many surfaces each node has seen so far
308*59599516SKenneth E. Jansen            do nn = 1, nshg
309*59599516SKenneth E. Jansen               if(nsurf(nn).eq.1)
310*59599516SKenneth E. Jansen     &              BCtmp(nn,4:6)=BCtmp(nn,4:6)+wnrm(nn,:)
311*59599516SKenneth E. Jansen               if(nsurf(nn).eq.2)
312*59599516SKenneth E. Jansen     &              BCtmp(nn,8:10)=BCtmp(nn,8:10)+wnrm(nn,:)
313*59599516SKenneth E. Jansen               if(nsurf(nn).gt.2)
314*59599516SKenneth E. Jansen     &              BCtmp(nn,4:6)=BCtmp(nn,4:6)+wnrm(nn,:)
315*59599516SKenneth E. Jansen            enddo
316*59599516SKenneth E. Jansenc That's all for this surface; move on to the next
317*59599516SKenneth E. Jansen         enddo                  ! end loop over surface ID's
318*59599516SKenneth E. Jansen
319*59599516SKenneth E. Jansenc
320*59599516SKenneth E. Jansenc.... complete the averaging, adjust iBC, adjust BCtmp
321*59599516SKenneth E. Jansenc
322*59599516SKenneth E. Jansen         wnrm(:,1:3)=zero
323*59599516SKenneth E. Jansen         do nn = 1, numnp       ! loop over all nodes
324*59599516SKenneth E. Jansenc leave nodes without wall-modeled surfaces unchanged
325*59599516SKenneth E. Jansen            if(nsurf(nn).eq.0) cycle
326*59599516SKenneth E. Jansenc
327*59599516SKenneth E. Jansenc mark this as a node with comp3
328*59599516SKenneth E. Jansenc
329*59599516SKenneth E. Jansen            ikp=0
330*59599516SKenneth E. Jansen            if(ibits(iBC(nn),3,3).eq.7 ) ikp=1
331*59599516SKenneth E. Jansenc         if(ibits(ibc(nn),3,3).eq.7 .and. BCtmp(nn,7).eq.zero) cycle
332*59599516SKenneth E. Jansenc find out which velocity BC's were set by the user, and cancel them
333*59599516SKenneth E. Jansen            ixset=ibits(iBC(nn),3,1)
334*59599516SKenneth E. Jansen            iyset=ibits(iBC(nn),4,1)
335*59599516SKenneth E. Jansen            izset=ibits(iBC(nn),5,1)
336*59599516SKenneth E. Jansen            iBC(nn)=iBC(nn)-ixset*8-iyset*16-izset*32
337*59599516SKenneth E. Jansenc
338*59599516SKenneth E. Jansenc save this stripped iBC for later un-fixing
339*59599516SKenneth E. Jansenc
340*59599516SKenneth E. Jansen            iBCSAV=iBC(nn)
341*59599516SKenneth E. Jansenc
342*59599516SKenneth E. Jansen            if(abs(itwmod).eq.1) then ! slip velocity wall-model
343*59599516SKenneth E. Jansenc If we're using the slip-velocity wall-model, then the velocity
344*59599516SKenneth E. Jansenc boundary condition for this node will depend on how many modeled
345*59599516SKenneth E. Jansenc surfaces share it...
346*59599516SKenneth E. Jansen               if(nsurf(nn).eq.1) then ! 1 modeled wall
347*59599516SKenneth E. Jansenc   If this node is part of only one modeled wall, then the component
348*59599516SKenneth E. Jansenc   of the velocity normal to the surface is set to zero.  In this case
349*59599516SKenneth E. Jansenc   only the first vector of BCtmp received normal contributions
350*59599516SKenneth E. Jansen                  iBC(nn)=iBC(nn)+8 ! assume normal is x-dominated first
351*59599516SKenneth E. Jansen                  wnrm(nn,:)=BCtmp(nn,4:6)
352*59599516SKenneth E. Jansen                  if(abs(wnrm(nn,3)).gt.abs(wnrm(nn,2)))then ! z beats y
353*59599516SKenneth E. Jansen                     if(abs(wnrm(nn,3)).gt.abs(wnrm(nn,1)))then ! z beats x
354*59599516SKenneth E. Jansen                        iBC(nn)=iBC(nn)+24
355*59599516SKenneth E. Jansen                     endif      ! z beats x
356*59599516SKenneth E. Jansen                  endif         ! z beats y
357*59599516SKenneth E. Jansen                  if(abs(wnrm(nn,2)).ge.abs(wnrm(nn,3)))then ! y beats z
358*59599516SKenneth E. Jansen                     if(abs(wnrm(nn,2)).gt.abs(wnrm(nn,1)))then ! y beats x
359*59599516SKenneth E. Jansen                        iBC(nn)=iBC(nn)+8
360*59599516SKenneth E. Jansen                     endif      ! y beats x
361*59599516SKenneth E. Jansen                  endif         ! y beats z           !(adjusted iBC)
362*59599516SKenneth E. Jansen                  BCtmp(nn,7)=zero
363*59599516SKenneth E. Jansen               else if(nsurf(nn).eq.2) then
364*59599516SKenneth E. Jansenc   If this node is shared by two modeled walls, then each wall
365*59599516SKenneth E. Jansenc   provides a normal vector along which the velocity must be zero.
366*59599516SKenneth E. Jansenc   This leaves only one "free" direction, along which the flow is nonzero.
367*59599516SKenneth E. Jansenc   The two normal vectors define a plane.  Any pair of non-parallel
368*59599516SKenneth E. Jansenc   vectors in this plane can also be specified, leaving the same free
369*59599516SKenneth E. Jansenc   direction.  Area-weighted average normal vectors for the two surfaces
370*59599516SKenneth E. Jansenc   sharing this node have been stored in BCtmp(4:6) and BCtmp(8:10)
371*59599516SKenneth E. Jansenc   We normalize the first of these, and then orthogonalize the second
372*59599516SKenneth E. Jansenc   against the first.  After then normalizing the second, we choose which
373*59599516SKenneth E. Jansenc   cartesian direction dominates each of the vectors, and adjust iBC
374*59599516SKenneth E. Jansenc   and BCtmp accordingly
375*59599516SKenneth E. Jansen                  e1=BCtmp(nn,4:6)
376*59599516SKenneth E. Jansen                  wmag=e1(1)*e1(1)+e1(2)*e1(2)+e1(3)*e1(3)
377*59599516SKenneth E. Jansen                  wmag=1./sqrt(wmag)
378*59599516SKenneth E. Jansen                  e1=e1*wmag    ! first vector is normalized
379*59599516SKenneth E. Jansenc
380*59599516SKenneth E. Jansen                  e2=BCtmp(nn,8:10)
381*59599516SKenneth E. Jansen                  wmag=e2(1)*e1(1)+e2(2)*e1(2)+e2(3)*e1(3) ! wmag=e2.n1
382*59599516SKenneth E. Jansen                  e2=e2-wmag*e1 ! second vector is orthogonalized
383*59599516SKenneth E. Jansenc
384*59599516SKenneth E. Jansen                  wmag=e2(1)*e2(1)+e2(2)*e2(2)+e2(3)*e2(3)
385*59599516SKenneth E. Jansen                  wmag=1./sqrt(wmag)
386*59599516SKenneth E. Jansen                  e2=e2*wmag    ! second vector is normalized
387*59599516SKenneth E. Jansenc
388*59599516SKenneth E. Jansen               ! idom tells us which component is currently dominant
389*59599516SKenneth E. Jansen               ! ivec(n) tells us which vector is dominated by comp n
390*59599516SKenneth E. Jansen                  ivec(:)=0     ! initialize
391*59599516SKenneth E. Jansen                  idom=1        ! assume x-comp dominates e1
392*59599516SKenneth E. Jansen                  bigcomp=abs(e1(1))
393*59599516SKenneth E. Jansen                  ivec(idom)=1
394*59599516SKenneth E. Jansen                  do i = 2, nsd
395*59599516SKenneth E. Jansen                     if(abs(e1(i)).gt.bigcomp) then
396*59599516SKenneth E. Jansen                        ivec(idom)=0
397*59599516SKenneth E. Jansen                        bigcomp=abs(e1(i))
398*59599516SKenneth E. Jansen                        idom=i
399*59599516SKenneth E. Jansen                        ivec(i)=1
400*59599516SKenneth E. Jansen                     endif
401*59599516SKenneth E. Jansen                  enddo
402*59599516SKenneth E. Jansen                  if(idom.ne.1) then
403*59599516SKenneth E. Jansen                     idom=1     ! assume x-comp dominates e2
404*59599516SKenneth E. Jansen                  else
405*59599516SKenneth E. Jansen                     idom=3     ! assume z-comp dominates e2
406*59599516SKenneth E. Jansen                  endif
407*59599516SKenneth E. Jansen                  bigcomp=abs(e2(idom))
408*59599516SKenneth E. Jansen
409*59599516SKenneth E. Jansen                  ivec(idom)=2
410*59599516SKenneth E. Jansen                  do i = 1, nsd
411*59599516SKenneth E. Jansen                     if(abs(e2(i)).gt.bigcomp) then
412*59599516SKenneth E. Jansen                        if(ivec(i).eq.1) cycle
413*59599516SKenneth E. Jansen                        ivec(idom)=0
414*59599516SKenneth E. Jansen                        bigcomp=abs(e2(i))
415*59599516SKenneth E. Jansen                        idom=i
416*59599516SKenneth E. Jansen                        ivec(i)=2
417*59599516SKenneth E. Jansen                     endif
418*59599516SKenneth E. Jansen                  enddo
419*59599516SKenneth E. Jansen               ! now we know which components dominate each vector
420*59599516SKenneth E. Jansen                  ixset=0       !
421*59599516SKenneth E. Jansen                  iyset=0       ! initialize
422*59599516SKenneth E. Jansen                  izset=0       !
423*59599516SKenneth E. Jansen                  if(ivec(1).ne.0) ixset=1
424*59599516SKenneth E. Jansen                  if(ivec(2).ne.0) iyset=1
425*59599516SKenneth E. Jansen                  if(ivec(3).ne.0) izset=1
426*59599516SKenneth E. Jansen                  ncomp=ixset+iyset+izset
427*59599516SKenneth E. Jansen                  if(ncomp.ne.2) write(*,*) 'WARNING: TROUBLE IN GENBC'
428*59599516SKenneth E. Jansen                  BCvecs(1,1:3)=e1(:)
429*59599516SKenneth E. Jansen                  BCvecs(2,1:3)=e2(:)
430*59599516SKenneth E. Jansen                  if((ixset.eq.1).and.(iyset.eq.1)) then
431*59599516SKenneth E. Jansen                     BCtmp(nn,4:6)=BCvecs(ivec(1),1:3)
432*59599516SKenneth E. Jansen                     BCtmp(nn,8:10)=BCvecs(ivec(2),1:3)
433*59599516SKenneth E. Jansen                  else if((ixset.eq.1).and.(izset.eq.1)) then
434*59599516SKenneth E. Jansen                     BCtmp(nn,4:6)=BCvecs(ivec(1),1:3)
435*59599516SKenneth E. Jansen                     BCtmp(nn,8:10)=BCvecs(ivec(3),1:3)
436*59599516SKenneth E. Jansen                  else if((iyset.eq.1).and.(izset.eq.1)) then
437*59599516SKenneth E. Jansen                     BCtmp(nn,4:6)=BCvecs(ivec(2),1:3)
438*59599516SKenneth E. Jansen                     BCtmp(nn,8:10)=BCvecs(ivec(3),1:3)
439*59599516SKenneth E. Jansen                  endif
440*59599516SKenneth E. Jansen                  iBC(nn)=iBC(nn)+ixset*8+iyset*16+izset*32
441*59599516SKenneth E. Jansen                  BCtmp(nn,7)=zero
442*59599516SKenneth E. Jansen                  BCtmp(nn,11)=zero
443*59599516SKenneth E. Jansen                  wnrm(nn,:)=BCtmp(nn,4:6)+BCtmp(nn,8:10)
444*59599516SKenneth E. Jansen               else if(nsurf(nn).gt.2) then
445*59599516SKenneth E. Jansenc     If this node is shared by more than two modeled walls, then
446*59599516SKenneth E. Jansenc     it is a corner node and it will be treated as having no slip
447*59599516SKenneth E. Jansenc     The normals for all surfaces beyond the first two were
448*59599516SKenneth E. Jansenc     contributed to the first vector of BCtmp
449*59599516SKenneth E. Jansen                  wnrm(nn,:)=BCtmp(nn,4:6)+BCtmp(nn,8:10)
450*59599516SKenneth E. Jansen                  iBC(nn)=iBC(nn)+56
451*59599516SKenneth E. Jansen                  BCtmp(nn,7)=zero
452*59599516SKenneth E. Jansen               endif            ! curved surface
453*59599516SKenneth E. Jansen            else if(abs(itwmod).eq.2) then ! effective viscosity wall-model
454*59599516SKenneth E. Jansenc Otherwise, we're using the effective-viscosity wall-model and the
455*59599516SKenneth E. Jansenc nodes on modeled surfaces are treated as no-slip
456*59599516SKenneth E. Jansen               iBC(nn)=iBC(nn)+56 ! set 3-comp
457*59599516SKenneth E. Jansen               if(itwmod.eq.-2) BCtmp(nn,7)=zero ! no slip
458*59599516SKenneth E. Jansen               if(nsurf(nn).eq.1)
459*59599516SKenneth E. Jansen     &              wnrm(nn,:)=BCtmp(nn,4:6)
460*59599516SKenneth E. Jansen               if(nsurf(nn).ge.2)
461*59599516SKenneth E. Jansen     &              wnrm(nn,:)=BCtmp(nn,4:6)+BCtmp(nn,8:10)
462*59599516SKenneth E. Jansen            endif
463*59599516SKenneth E. Jansenc Now normalize the wall normal for this node
464*59599516SKenneth E. Jansen            if(itwmod.ne.0) then
465*59599516SKenneth E. Jansen               wmag=sqrt(wnrm(nn,1)*wnrm(nn,1)
466*59599516SKenneth E. Jansen     &              +wnrm(nn,2)*wnrm(nn,2)+wnrm(nn,3)*wnrm(nn,3))
467*59599516SKenneth E. Jansen               wnrm(nn,:)=wnrm(nn,:)/wmag
468*59599516SKenneth E. Jansen            endif
469*59599516SKenneth E. Jansenc
470*59599516SKenneth E. Jansenc.... put back the comp3 info for bctmp
471*59599516SKenneth E. Jansenc
472*59599516SKenneth E. Jansen            if(ikp.eq.1) then
473*59599516SKenneth E. Jansen               iBC(nn)=iBCSAV+56 ! put it back to a comp3
474*59599516SKenneth E. Jansen               BCtmp(nn,:)=BCtmpSAV(nn,:) ! ditto
475*59599516SKenneth E. Jansen            endif
476*59599516SKenneth E. Jansen         enddo                  ! loop over all nodes
477*59599516SKenneth E. Jansen      endif                     ! end "if there are any surfID's"
478*59599516SKenneth E. Jansenc
479*59599516SKenneth E. Jansenc  If you are using the turbulence wall with axisymmetry you
480*59599516SKenneth E. Jansenc  need to modify the axisymmetry angle to account for the discrete
481*59599516SKenneth E. Jansenc  normals at the wall being different than the exact normals
482*59599516SKenneth E. Jansenc
483*59599516SKenneth E. Jansenc find the my normal, my partners normal and correct the angle
484*59599516SKenneth E. Jansenc$$$
485*59599516SKenneth E. Jansenc$$$        do i=1,numnp
486*59599516SKenneth E. Jansenc$$$           wmag = wnrm(i,1) * wnrm(i,1)
487*59599516SKenneth E. Jansenc$$$     &          + wnrm(i,2) * wnrm(i,2)
488*59599516SKenneth E. Jansenc$$$     &          + wnrm(i,3) * wnrm(i,3)
489*59599516SKenneth E. Jansenc$$$c
490*59599516SKenneth E. Jansenc$$$c  only "fix" wall nodes...other nodes still have a zero in wnrm
491*59599516SKenneth E. Jansenc$$$c
492*59599516SKenneth E. Jansenc$$$           if((btest(iBC(i),12)).and.(wmag.ne.zero)) then
493*59599516SKenneth E. Jansenc$$$              BCtmp(i,1)=acos( wnrm(i,1) * wnrm(iper(i),1)
494*59599516SKenneth E. Jansenc$$$     &                       + wnrm(i,2) * wnrm(iper(i),2)
495*59599516SKenneth E. Jansenc$$$     &                       + wnrm(i,3) * wnrm(iper(i),3) )
496*59599516SKenneth E. Jansenc$$$           endif
497*59599516SKenneth E. Jansenc$$$        enddo
498*59599516SKenneth E. Jansenc
499*59599516SKenneth E. Jansenc.... return
500*59599516SKenneth E. Jansenc
501*59599516SKenneth E. Jansen      return
502*59599516SKenneth E. Jansenc
503*59599516SKenneth E. Jansen      end
504*59599516SKenneth E. Jansen
505*59599516SKenneth E. Jansen
506*59599516SKenneth E. Jansen      subroutine gensidcount(nsidg)
507*59599516SKenneth E. Jansenc---------------------------------------------------------------------
508*59599516SKenneth E. Jansenc
509*59599516SKenneth E. Jansenc This routine counts up the total number of surface ID's across
510*59599516SKenneth E. Jansenc all processors and makes a list of them
511*59599516SKenneth E. Jansenc
512*59599516SKenneth E. Jansenc Inputs:
513*59599516SKenneth E. Jansenc     iBCB        natural boundary condition switches and surfIDs
514*59599516SKenneth E. Jansenc
515*59599516SKenneth E. Jansenc Outputs:
516*59599516SKenneth E. Jansenc     nsidg       number of surface ID's globally (including all procs)
517*59599516SKenneth E. Jansenc     sidmapg     global list of surface ID's, lowest to highest
518*59599516SKenneth E. Jansenc
519*59599516SKenneth E. Jansenc---------------------------------------------------------------------
520*59599516SKenneth E. Jansen      use pointer_data ! access to miBCB
521*59599516SKenneth E. Jansen      use turbSA ! access to sidmapg
522*59599516SKenneth E. Jansen      include "common.h"
523*59599516SKenneth E. Jansen      include "mpif.h"
524*59599516SKenneth E. Jansenc
525*59599516SKenneth E. Jansen      integer newflag, i, j
526*59599516SKenneth E. Jansen      integer nsidl             ! number of surface IDs on-proc
527*59599516SKenneth E. Jansen      integer nsidt             ! sum of all nsidl's
528*59599516SKenneth E. Jansen      integer nsidg             ! number of surface IDs globally
529*59599516SKenneth E. Jansen      integer nsid(numpe)       ! array of nsidl's
530*59599516SKenneth E. Jansen      integer idisp(numpe)      ! needed by mpi: see note below
531*59599516SKenneth E. Jansen      type llnod                ! structure for one node in a linked list
532*59599516SKenneth E. Jansen        integer :: value
533*59599516SKenneth E. Jansen        type (llnod), pointer :: next
534*59599516SKenneth E. Jansen      end type llnod
535*59599516SKenneth E. Jansen      type (llnod), pointer :: sidlist ! points to first elt of linked list
536*59599516SKenneth E. Jansen      type (llnod), pointer :: sidelt  ! points to generic elt of linked list
537*59599516SKenneth E. Jansen      integer, allocatable :: sidmapl(:) ! list of surfID's on-proc
538*59599516SKenneth E. Jansen      integer, allocatable :: temp(:)    ! temp space
539*59599516SKenneth E. Jansenc      integer iBCB(numelb,ndiBCB) ! column 1: naturalBC switches
540*59599516SKenneth E. Jansen                                   ! column 2: surface ID's
541*59599516SKenneth E. Jansenc Since we don't know a priori how many surface ID's there are,
542*59599516SKenneth E. Jansenc on-processor or globally, we will store the ID's as an open-ended
543*59599516SKenneth E. Jansenc link list while we determine the total number of distinct ID's
544*59599516SKenneth E. Jansen      allocate (sidlist) ! allocate the first element of the sid
545*59599516SKenneth E. Jansen                         ! linked list and point sidlist to it
546*59599516SKenneth E. Jansen      nsidl=0            ! initialize # of sids on this processor
547*59599516SKenneth E. Jansen      nsidg=0
548*59599516SKenneth E. Jansen      do iblk=1, nelblb  ! loop over boundary element blocks
549*59599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
550*59599516SKenneth E. Jansen         do i = 1, npro         ! loop over boundary elements (belts)
551*59599516SKenneth E. Jansen            iBCB2=miBCB(iblk)%p(i,2)
552*59599516SKenneth E. Jansen            if(iBCB2.ne.zero) then ! if a sid is given for this belt
553*59599516SKenneth E. Jansen               if(nsidl.eq.0) then     !   if this is the first sid we've seen
554*59599516SKenneth E. Jansen                  nsidl=1              !       increment our count and
555*59599516SKenneth E. Jansen                  sidlist % value = iBCB2    ! record its value
556*59599516SKenneth E. Jansen               else                    !   if this isn't the first sid
557*59599516SKenneth E. Jansen                  newflag = 1          !     assume this is a new sid
558*59599516SKenneth E. Jansen                  sidelt => sidlist    !     starting with the first sid
559*59599516SKenneth E. Jansen                  do j = 1, nsidl      !     check the assumption
560*59599516SKenneth E. Jansen                     if((sidelt % value).eq.iBCB2) then
561*59599516SKenneth E. Jansen                        newflag=0      !        ...
562*59599516SKenneth E. Jansen                     endif
563*59599516SKenneth E. Jansen                     if(j.ne.nsidl) then !      ...
564*59599516SKenneth E. Jansen                        sidelt => sidelt % next
565*59599516SKenneth E. Jansen                     endif!                     ...
566*59599516SKenneth E. Jansen                  enddo
567*59599516SKenneth E. Jansen                  if(newflag.eq.1) then!     if it really is new to us
568*59599516SKenneth E. Jansen                     nsidl = nsidl + 1 !         increment our count
569*59599516SKenneth E. Jansen                     allocate (sidelt % next)!   tack a new element to the end
570*59599516SKenneth E. Jansen                     sidelt => sidelt % next!    point to the new element
571*59599516SKenneth E. Jansen                     sidelt % value = iBCB2    ! record the new sid
572*59599516SKenneth E. Jansen                  endif ! (really is a new sid)
573*59599516SKenneth E. Jansen               endif ! (first sid)
574*59599516SKenneth E. Jansen            endif ! (belt has a sid)
575*59599516SKenneth E. Jansen         enddo ! (loop over belts)
576*59599516SKenneth E. Jansen      enddo ! (loop over boundary element blocks)
577*59599516SKenneth E. Jansenc Copy the data from the linked list to a more easily-used array form
578*59599516SKenneth E. Jansen      if(nsidl.gt.0) then
579*59599516SKenneth E. Jansen         allocate( sidmapl(nsidl) )
580*59599516SKenneth E. Jansen         sidelt => sidlist      !     starting with the first sid
581*59599516SKenneth E. Jansen         do j = 1, nsidl
582*59599516SKenneth E. Jansen            sidmapl(j)=sidelt%value
583*59599516SKenneth E. Jansen            if(j.ne.nsidl) sidelt => sidelt%next
584*59599516SKenneth E. Jansen         enddo
585*59599516SKenneth E. Jansen      else
586*59599516SKenneth E. Jansen	allocate( sidmapl(1)) ! some compilers/MPI don't like to send unallocated arrays
587*59599516SKenneth E. Jansen      endif
588*59599516SKenneth E. Jansenc Gather the number of surface ID's on each processor
589*59599516SKenneth E. Jansen      if (numpe.gt.1) then      ! multiple processors
590*59599516SKenneth E. Jansenc write the number of surfID's on the jth processor to slot j of nsid
591*59599516SKenneth E. Jansen         call MPI_ALLGATHER(nsidl,1,MPI_INTEGER,nsid,1,
592*59599516SKenneth E. Jansen     &          MPI_INTEGER,MPI_COMM_WORLD,ierr)
593*59599516SKenneth E. Jansenc count up the total number of surfID's among all processes
594*59599516SKenneth E. Jansen         nsidt=0
595*59599516SKenneth E. Jansen         do j=1,numpe
596*59599516SKenneth E. Jansen            nsidt=nsidt+nsid(j)
597*59599516SKenneth E. Jansen         enddo
598*59599516SKenneth E. Jansen      else                      ! single processor
599*59599516SKenneth E. Jansenc the local information is the global information for single-processor
600*59599516SKenneth E. Jansen         nsid=nsidl
601*59599516SKenneth E. Jansen         nsidt=nsidl
602*59599516SKenneth E. Jansen      endif                     ! if-else for multiple processors
603*59599516SKenneth E. Jansen      if(nsidt.gt.0) then
604*59599516SKenneth E. Jansenc
605*59599516SKenneth E. Jansenc  Make all-processor surfID collage
606*59599516SKenneth E. Jansenc
607*59599516SKenneth E. Jansenc there will be some duplicate surface ID's when we gather, so we
608*59599516SKenneth E. Jansenc will use a temporary array
609*59599516SKenneth E. Jansen         allocate (temp(nsidt))
610*59599516SKenneth E. Jansen         if (numpe.gt.1) then   ! multiple processors
611*59599516SKenneth E. Jansenc we will gather surfID's from local on-proc sets to a global set
612*59599516SKenneth E. Jansenc we will stack each processor's surfID list atop that of the previous
613*59599516SKenneth E. Jansenc processor.  If the list for processor i is called sidmapi, then our
614*59599516SKenneth E. Jansenc global coordinate list sidmap will look like this:
615*59599516SKenneth E. Jansenc ---------------------------------------------------------------
616*59599516SKenneth E. Jansenc | sidmap1       | sidmap2           | sidmap3   |   ...       |
617*59599516SKenneth E. Jansenc ---------------------------------------------------------------
618*59599516SKenneth E. Jansenc  <---nsid(1)---> <-----nsid(2)-----> <-nsid(3)->
619*59599516SKenneth E. Jansenc  <------------------------nsidt-----------------------...---->
620*59599516SKenneth E. Jansenc To accomplish this with MPI, we use MPI_ALLGATHERV, summarized as:
621*59599516SKenneth E. JansencMPI_ALLGATHERV(sendbuf,sendcount,sendtype,recvbuf,recvcount,disp,recvtype,comm)
622*59599516SKenneth E. Jansenc[ IN sendbuf] starting address of send buffer (choice)
623*59599516SKenneth E. Jansenc[ IN sendcount] number of elements in send buffer (integer)
624*59599516SKenneth E. Jansenc[ IN sendtype] data type of send buffer elements (handle)
625*59599516SKenneth E. Jansenc[ OUT recvbuf] address of receive buffer (choice)
626*59599516SKenneth E. Jansenc[ IN recvcount] number of elements received from each process (int array)
627*59599516SKenneth E. Jansenc[ IN disp] displacement array
628*59599516SKenneth E. Jansenc[ IN recvtype] data type of receive buffer elements (handle)
629*59599516SKenneth E. Jansenc[ IN comm] communicator (handle)
630*59599516SKenneth E. Jansenc The displacement array disp is an array of integers in which the jth
631*59599516SKenneth E. Jansenc entry indicates which slot of sidmap marks beginning of sidmapj
632*59599516SKenneth E. Jansenc So, first we will build this displacement array
633*59599516SKenneth E. Jansen            idisp(:)=0      ! starting with zero, since MPI likes C-numbering
634*59599516SKenneth E. Jansen            do j=2,numpe
635*59599516SKenneth E. Jansen               idisp(j)=idisp(j-1)+nsid(j-1) ! see diagram above
636*59599516SKenneth E. Jansen            enddo
637*59599516SKenneth E. Jansenc Now, we gather the data
638*59599516SKenneth E. Jansen            call MPI_ALLGATHERV(sidmapl(:),nsidl,
639*59599516SKenneth E. Jansen     &           MPI_INTEGER,temp(:),nsid,idisp,
640*59599516SKenneth E. Jansen     &           MPI_INTEGER,MPI_COMM_WORLD,ierr)
641*59599516SKenneth E. Jansenc sort surfID's, lowest to highest
642*59599516SKenneth E. Jansen            isorted = 0
643*59599516SKenneth E. Jansen            do while (isorted.eq.0) ! while the list isn't sorted
644*59599516SKenneth E. Jansen               isorted = 1      ! assume the list is sorted this time
645*59599516SKenneth E. Jansen               do j = 2, nsidt  ! loop over ID's
646*59599516SKenneth E. Jansen                  if(temp(j).lt.temp(j-1)) then ! ID exceeds predecessor
647*59599516SKenneth E. Jansen                     itmp=temp(j-1)
648*59599516SKenneth E. Jansen                     temp(j-1)=temp(j)
649*59599516SKenneth E. Jansen                     temp(j)=itmp
650*59599516SKenneth E. Jansen                     isorted=0  ! not sorted yet
651*59599516SKenneth E. Jansen                  endif
652*59599516SKenneth E. Jansen               enddo            !loop over ID's
653*59599516SKenneth E. Jansen            enddo               ! while not sorted
654*59599516SKenneth E. Jansenc determine the total number of distinct surfID's, globally
655*59599516SKenneth E. Jansen            nsidg=nsidt         ! assume there are no duplicate SurfID's
656*59599516SKenneth E. Jansen            do j = 2, nsidt
657*59599516SKenneth E. Jansen               if(temp(j).eq.temp(j-1)) nsidg = nsidg - 1 ! correction
658*59599516SKenneth E. Jansen            enddo
659*59599516SKenneth E. Jansenc create duplicate-free surfID list
660*59599516SKenneth E. Jansen            allocate( sidmapg(nsidg) )
661*59599516SKenneth E. Jansen            sidmapg(1)=temp(1)
662*59599516SKenneth E. Jansen            nsidg = 1
663*59599516SKenneth E. Jansen            do j = 2, nsidt
664*59599516SKenneth E. Jansen               if(temp(j).ne.temp(j-1)) then
665*59599516SKenneth E. Jansen                  nsidg = nsidg + 1
666*59599516SKenneth E. Jansen                  sidmapg(nsidg)=temp(j)
667*59599516SKenneth E. Jansen               endif
668*59599516SKenneth E. Jansen            enddo
669*59599516SKenneth E. Jansen            deallocate( temp )
670*59599516SKenneth E. Jansen         else                   ! single-processor
671*59599516SKenneth E. Jansenc global data is local data in single processor case
672*59599516SKenneth E. Jansen            nsidg=nsidl
673*59599516SKenneth E. Jansen            allocate( sidmapg(nsidg) )
674*59599516SKenneth E. Jansen            sidmapg=sidmapl
675*59599516SKenneth E. Jansen         endif                  ! if-else multiple processor
676*59599516SKenneth E. Jansenc
677*59599516SKenneth E. Jansen      endif                     ! if at least one surfid
678*59599516SKenneth E. Jansenc
679*59599516SKenneth E. Jansen      return
680*59599516SKenneth E. Jansenc
681*59599516SKenneth E. Jansen      end
682*59599516SKenneth E. Jansen
683*59599516SKenneth E. Jansen
684*59599516SKenneth E. Jansen
685*59599516SKenneth E. Jansen      subroutine genotwn(x, BCtmp, iBC, nsurf)
686*59599516SKenneth E. Jansenc
687*59599516SKenneth E. Jansenc----------------------------------------------------------------------
688*59599516SKenneth E. Jansenc  This routine determines which node to use as the first node off the
689*59599516SKenneth E. Jansenc  wall in near-wall modeling traction calculations for each wall node.
690*59599516SKenneth E. Jansenc  Each wall node has a normal, calculated from the wall elements to
691*59599516SKenneth E. Jansenc  which that node belongs.  We seek the node within the boundary
692*59599516SKenneth E. Jansenc  element that lies closest to the line defined by the normal vector.
693*59599516SKenneth E. Jansenc  We create normalized vectors pointing from the wall node in
694*59599516SKenneth E. Jansenc  question to each of the nodes in the boundary element. The vector
695*59599516SKenneth E. Jansenc  that has the largest projection onto the wall node's normal vector
696*59599516SKenneth E. Jansenc  points to the node we want.  Nodes that are not on a wall point to
697*59599516SKenneth E. Jansenc  themselves as their own off-the-wall node.
698*59599516SKenneth E. Jansenc
699*59599516SKenneth E. Jansenc input:
700*59599516SKenneth E. Jansenc  x     (nshg,3)     :        : nodal position vectors
701*59599516SKenneth E. Jansenc  wnrm  (nshg,3)     : (mod)  : normal vectors for each node
702*59599516SKenneth E. Jansenc  iBCB5 (nshg)       : (file) : wall flag for belt
703*59599516SKenneth E. Jansenc  ienb  (numelb,nen) : (file) : connectivity for belts
704*59599516SKenneth E. Jansenc
705*59599516SKenneth E. Jansenc output:
706*59599516SKenneth E. Jansenc  otwn  (nshg)       : (mod)  : off-the-wall nodes for each node
707*59599516SKenneth E. Jansenc
708*59599516SKenneth E. Jansenc
709*59599516SKenneth E. Jansenc Kenneth Jansen, Summer 2000 (algorithm)
710*59599516SKenneth E. Jansenc Michael Yaworski, Summer 2000 (code)
711*59599516SKenneth E. Jansenc----------------------------------------------------------------------
712*59599516SKenneth E. Jansenc
713*59599516SKenneth E. Jansen      use pointer_data          ! used for mienb, miBCB
714*59599516SKenneth E. Jansen      use turbSA
715*59599516SKenneth E. Jansen      include "common.h"
716*59599516SKenneth E. Jansenc
717*59599516SKenneth E. Jansen      integer iel, nod, can
718*59599516SKenneth E. Jansen      real*8 vec(3), leng, dp, bigdp, lil
719*59599516SKenneth E. Jansen      real*8 x(numnp,nsd),BCtmp(nshg,ndof+7)
720*59599516SKenneth E. Jansen      integer iBC(nshg), nsurf(nshg)
721*59599516SKenneth E. Jansen      integer gbits
722*59599516SKenneth E. Jansen      integer, allocatable :: ienb(:)
723*59599516SKenneth E. Jansen
724*59599516SKenneth E. Jansen      allocate( otwn(nshg) )
725*59599516SKenneth E. Jansenc
726*59599516SKenneth E. Jansenc initialize otwn to oneself
727*59599516SKenneth E. Jansenc
728*59599516SKenneth E. Jansen      do nod = 1, nshg
729*59599516SKenneth E. Jansen         otwn(nod)=nod
730*59599516SKenneth E. Jansen      enddo
731*59599516SKenneth E. Jansenc
732*59599516SKenneth E. Jansenc determine otwn for each wall node
733*59599516SKenneth E. Jansenc
734*59599516SKenneth E. Jansen      do iblk=1, nelblb         ! loop over boundary element blocks
735*59599516SKenneth E. Jansen         npro = lcblkb(1,iblk+1)-lcblkb(1,iblk)
736*59599516SKenneth E. Jansen         nenl  = lcblkb(5,iblk)
737*59599516SKenneth E. Jansen         nenbl = lcblkb(6,iblk)
738*59599516SKenneth E. Jansen         nshl = lcblkb(9,iblk)
739*59599516SKenneth E. Jansen         allocate( ienb(nshl) )
740*59599516SKenneth E. Jansen         do i = 1, npro         ! loop over boundary elements
741*59599516SKenneth E. Jansen            iBCB1=miBCB(iblk)%p(i,1)
742*59599516SKenneth E. Jansen            iBCB2=miBCB(iblk)%p(i,2)
743*59599516SKenneth E. Jansen            ienb(1:nshl)=mienb(iblk)%p(i,1:nshl)
744*59599516SKenneth E. Jansen            if (btest(iBCB1,4)) then ! (on the wall)
745*59599516SKenneth E. Jansen               do nod = 1, nenbl !   for each wall node in this belt
746*59599516SKenneth E. Jansen                  bigdp = zero  !     initialize largest dot product
747*59599516SKenneth E. Jansen                  do can=nenbl+1,nenl !  loop over off-the-wall candidates
748*59599516SKenneth E. Jansen                     nn=ienb(can)
749*59599516SKenneth E. Jansen                               !       don't bother with wall nodes
750*59599516SKenneth E. Jansen                     if(nsurf(nn).ne.0) cycle
751*59599516SKenneth E. Jansen                               !       don't bother with no-slip nodes
752*59599516SKenneth E. Jansen                     if(ibits(iBC(nn),3,3).eq.7 .and.
753*59599516SKenneth E. Jansen     &                    BCtmp(nn,7).eq.zero) cycle
754*59599516SKenneth E. Jansenc                              !       calculate candidate vector
755*59599516SKenneth E. Jansen                     vec(:)=x(ienb(can),:)-x(ienb(nod),:)
756*59599516SKenneth E. Jansen                     leng=sqrt(vec(1)**2+vec(2)**2+vec(3)**2)
757*59599516SKenneth E. Jansen                     vec(:)=vec(:)/leng
758*59599516SKenneth E. Jansenc                              !       calculate dot product with wnrm
759*59599516SKenneth E. Jansen                     vec(:)=vec(:)*wnrm(ienb(nod),:)
760*59599516SKenneth E. Jansen                     dp=vec(1)+vec(2)+vec(3)
761*59599516SKenneth E. Jansenc                              !       set candidate as otwn if necessary
762*59599516SKenneth E. Jansenc                              !       (wnrm points into fluid, so
763*59599516SKenneth E. Jansenc                              !        we want the most positive dp)
764*59599516SKenneth E. Jansen                     if (dp.gt.bigdp) then
765*59599516SKenneth E. Jansen                        otwn(ienb(nod)) = ienb(can)
766*59599516SKenneth E. Jansen                        bigdp=dp
767*59599516SKenneth E. Jansen                     endif
768*59599516SKenneth E. Jansen                  enddo         !(loop over off-the-wall candidates)
769*59599516SKenneth E. Jansen               enddo            !(loop over wall nodes in current belt)
770*59599516SKenneth E. Jansen            endif
771*59599516SKenneth E. Jansen         enddo                  !(loop over elts in block)
772*59599516SKenneth E. Jansen         deallocate(ienb)
773*59599516SKenneth E. Jansen      enddo                     !(loop over boundary element blocks)
774*59599516SKenneth E. Jansen      do nn = 1, nshg
775*59599516SKenneth E. Jansen         if((otwn(nn).eq.nn).and.(nsurf(nn).ne.0)) then ! if a node on a
776*59599516SKenneth E. Jansen                                                    ! modeled surface
777*59599516SKenneth E. Jansen                                                    ! didn't find a node
778*59599516SKenneth E. Jansen                                                    ! off the wall
779*59599516SKenneth E. Jansen            lil = 1.0e32 ! distance to current closest prospect
780*59599516SKenneth E. Jansen            do can = 1, nshg ! loop over nodes
781*59599516SKenneth E. Jansen               if(nsurf(can).eq.0) then ! if this candidate is off the
782*59599516SKenneth E. Jansen                                        ! wall
783*59599516SKenneth E. Jansen                  if(ibits(iBC(can),3,3).eq.7 .and.
784*59599516SKenneth E. Jansen     &               BCtmp(can,7).eq.zero) then  ! no slip nodes not allowed
785*59599516SKenneth E. Jansen                  else          ! not a forbidden node
786*59599516SKenneth E. Jansen                     vec(:)=x(nn,:)-x(can,:)
787*59599516SKenneth E. Jansen                     leng=vec(1)**2+vec(2)**2+vec(3)**2
788*59599516SKenneth E. Jansen                     if(leng.lt.lil) then ! this candidate is closest so far
789*59599516SKenneth E. Jansen                        lil=leng
790*59599516SKenneth E. Jansen                        otwn(nn)=can
791*59599516SKenneth E. Jansen                     endif      ! closest so far
792*59599516SKenneth E. Jansen                  endif  ! end of no slip nodes
793*59599516SKenneth E. Jansen               endif ! off-the-wall check
794*59599516SKenneth E. Jansen            enddo ! loop over nodes
795*59599516SKenneth E. Jansen         endif ! lonely wall-node check
796*59599516SKenneth E. Jansen      enddo ! loop over nodes
797*59599516SKenneth E. Jansenc
798*59599516SKenneth E. Jansen      return
799*59599516SKenneth E. Jansenc
800*59599516SKenneth E. Jansen      end
801*59599516SKenneth E. Jansen
802