xref: /phasta/phSolver/common/asbnabi.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen      subroutine AsBNABI ( x,       shpb,
2*59599516SKenneth E. Jansen     &                   ienb,  iBCB)
3*59599516SKenneth E. Jansenc
4*59599516SKenneth E. Jansenc----------------------------------------------------------------------
5*59599516SKenneth E. Jansenc
6*59599516SKenneth E. Jansenc This routine computes and assembles the data corresponding to the
7*59599516SKenneth E. Jansenc  boundary elements.
8*59599516SKenneth E. Jansenc
9*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
10*59599516SKenneth E. Jansenc----------------------------------------------------------------------
11*59599516SKenneth E. Jansenc
12*59599516SKenneth E. Jansen      use pvsQbi
13*59599516SKenneth E. Jansen        include "common.h"
14*59599516SKenneth E. Jansenc
15*59599516SKenneth E. Jansen        dimension xlb(npro,nenl,nsd),    bnorm(npro,nsd),
16*59599516SKenneth E. Jansen     &            rl(npro,nshl,nsd),     WdetJb(npro),
17*59599516SKenneth E. Jansen     &            Wfactor(npro)
18*59599516SKenneth E. Jansen
19*59599516SKenneth E. Jansen        dimension x(numnp,nsd),
20*59599516SKenneth E. Jansen     &            shpb(nshl,ngaussb),      shglb(nsd,nshl,ngaussb),
21*59599516SKenneth E. Jansen     &            ienb(npro,nshl),
22*59599516SKenneth E. Jansen     &            iBCB(npro,ndiBCB)
23*59599516SKenneth E. Jansen
24*59599516SKenneth E. Jansen        dimension lnode(27),               sgn(npro,nshl),
25*59599516SKenneth E. Jansen     &            shpfun(npro,nshl),        shdrv(npro,nsd,nshl)
26*59599516SKenneth E. Jansen
27*59599516SKenneth E. Jansenc
28*59599516SKenneth E. Jansen        dimension dxdxib(npro,nsd,nsd),      temp(npro),
29*59599516SKenneth E. Jansen     &            temp1(npro),               temp2(npro),
30*59599516SKenneth E. Jansen     &            temp3(npro),
31*59599516SKenneth E. Jansen     &            v1(npro,nsd),              v2(npro,nsd)
32*59599516SKenneth E. Jansen
33*59599516SKenneth E. Jansenc
34*59599516SKenneth E. Jansenc.... get the matrix of mode signs for the hierarchic basis functions
35*59599516SKenneth E. Jansenc
36*59599516SKenneth E. Jansen        if (ipord .gt. 1) then
37*59599516SKenneth E. Jansen           call getsgn(ienb,sgn)
38*59599516SKenneth E. Jansen        endif
39*59599516SKenneth E. Jansenc
40*59599516SKenneth E. Jansenc.... gather the variables
41*59599516SKenneth E. Jansenc
42*59599516SKenneth E. Jansen        call localx(x,      xlb,    ienb,   nsd,    'gather  ')
43*59599516SKenneth E. Jansenc
44*59599516SKenneth E. Jansenc.... get the boundary element residuals
45*59599516SKenneth E. Jansenc
46*59599516SKenneth E. Jansen        rl  = zero
47*59599516SKenneth E. Jansenc
48*59599516SKenneth E. Jansenc.... compute the nodes which lie on the boundary (hierarchic)
49*59599516SKenneth E. Jansenc
50*59599516SKenneth E. Jansen        call getbnodes(lnode)
51*59599516SKenneth E. Jansenc
52*59599516SKenneth E. Jansenc.... compute the normal to the boundary. This is achieved by taking
53*59599516SKenneth E. Jansenc     the cross product of two vectors in the plane of the 2-d
54*59599516SKenneth E. Jansenc     boundary face.
55*59599516SKenneth E. Jansenc
56*59599516SKenneth E. Jansen  	if(lcsyst.eq.1) then	  ! set to curl into element all others out
57*59599516SKenneth E. Jansen  	   ipt2=2
58*59599516SKenneth E. Jansen  	   ipt3=3
59*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.2) then
60*59599516SKenneth E. Jansen  	   ipt2=4
61*59599516SKenneth E. Jansen  	   ipt3=2
62*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.3) then
63*59599516SKenneth E. Jansen  	   ipt2=3
64*59599516SKenneth E. Jansen  	   ipt3=2
65*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.4) then
66*59599516SKenneth E. Jansen  	   ipt2=2
67*59599516SKenneth E. Jansen  	   ipt3=4
68*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.5) then
69*59599516SKenneth E. Jansen  	   ipt2=4
70*59599516SKenneth E. Jansen  	   ipt3=2
71*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.6) then
72*59599516SKenneth E. Jansen  	   ipt2=2
73*59599516SKenneth E. Jansen  	   ipt3=5
74*59599516SKenneth E. Jansen  	endif
75*59599516SKenneth E. Jansen  	v1 = xlb(:,ipt2,:) - xlb(:,1,:)
76*59599516SKenneth E. Jansen  	v2 = xlb(:,ipt3,:) - xlb(:,1,:)
77*59599516SKenneth E. Jansenc
78*59599516SKenneth E. Jansenc compute cross product
79*59599516SKenneth E. Jansenc
80*59599516SKenneth E. Jansen        temp1 = v1(:,2) * v2(:,3) - v2(:,2) * v1(:,3)
81*59599516SKenneth E. Jansen        temp2 = v2(:,1) * v1(:,3) - v1(:,1) * v2(:,3)
82*59599516SKenneth E. Jansen        temp3 = v1(:,1) * v2(:,2) - v2(:,1) * v1(:,2)
83*59599516SKenneth E. Jansenc
84*59599516SKenneth E. Jansenc mag is area for quads, twice area for tris
85*59599516SKenneth E. Jansenc
86*59599516SKenneth E. Jansen  	temp	   = one / sqrt ( temp1**2 + temp2**2 + temp3**2 )
87*59599516SKenneth E. Jansen  	bnorm(:,1) = temp1 * temp
88*59599516SKenneth E. Jansen  	bnorm(:,2) = temp2 * temp
89*59599516SKenneth E. Jansen  	bnorm(:,3) = temp3 * temp
90*59599516SKenneth E. Jansenc
91*59599516SKenneth E. Jansen  	if (lcsyst .eq. 1) then
92*59599516SKenneth E. Jansen  	   Wfactor(:) = one / (four*temp(:))
93*59599516SKenneth E. Jansen  	elseif (lcsyst .eq. 2) then
94*59599516SKenneth E. Jansen  	   Wfactor(:) = one / (four*temp(:))
95*59599516SKenneth E. Jansen  	elseif (lcsyst .eq. 3) then
96*59599516SKenneth E. Jansen  	   Wfactor(:) = one / (two*temp(:))
97*59599516SKenneth E. Jansen  	elseif (lcsyst .eq. 4) then
98*59599516SKenneth E. Jansen  	   Wfactor(:) = one / (four*temp(:))
99*59599516SKenneth E. Jansen  	elseif (lcsyst .eq. 5) then
100*59599516SKenneth E. Jansen  	   Wfactor(:) = one / (four*temp(:))
101*59599516SKenneth E. Jansen  	elseif (lcsyst .eq. 6) then
102*59599516SKenneth E. Jansen  	   Wfactor(:) = one / (two*temp(:))
103*59599516SKenneth E. Jansen  	endif
104*59599516SKenneth E. Jansenc
105*59599516SKenneth E. Jansen        do iel=1,npro
106*59599516SKenneth E. Jansen           if (.not.btest(iBCB(iel,1),1)) then
107*59599516SKenneth E. Jansen              Wfactor(iel) = zero  ! we want zeros where we are not integrating
108*59599516SKenneth E. Jansen           endif
109*59599516SKenneth E. Jansen        enddo
110*59599516SKenneth E. Jansenc
111*59599516SKenneth E. Jansen        ngaussb = nintb(lcsyst)
112*59599516SKenneth E. Jansenc
113*59599516SKenneth E. Jansenc.... loop through the integration points
114*59599516SKenneth E. Jansenc
115*59599516SKenneth E. Jansen        do intp = 1, ngaussb
116*59599516SKenneth E. Jansenc
117*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point
118*59599516SKenneth E. Jansenc
119*59599516SKenneth E. Jansen           shglb=zero  ! protect debugger
120*59599516SKenneth E. Jansen           call getshpb(shpb,        shglb,        sgn,
121*59599516SKenneth E. Jansen     &              shpfun,       shdrv)
122*59599516SKenneth E. Jansenc
123*59599516SKenneth E. Jansen           WdetJb(:) = Qwtb(lcsyst,intp) * Wfactor(:)
124*59599516SKenneth E. Jansenc
125*59599516SKenneth E. Jansenc  Now lets calculate Integral N_(a:e)^i n_i d Gamma
126*59599516SKenneth E. Jansenc
127*59599516SKenneth E. Jansen           do n = 1, nshlb
128*59599516SKenneth E. Jansen              nodlcl = lnode(n)
129*59599516SKenneth E. Jansen              rl(:,nodlcl,1) = rl(:,nodlcl,1) +
130*59599516SKenneth E. Jansen     &             shpfun(:,nodlcl) * bnorm(:,1) * WdetJb(:)
131*59599516SKenneth E. Jansen              rl(:,nodlcl,2) = rl(:,nodlcl,2) +
132*59599516SKenneth E. Jansen     &             shpfun(:,nodlcl) * bnorm(:,2) * WdetJb(:)
133*59599516SKenneth E. Jansen              rl(:,nodlcl,3) = rl(:,nodlcl,3) +
134*59599516SKenneth E. Jansen     &             shpfun(:,nodlcl) * bnorm(:,3) * WdetJb(:)
135*59599516SKenneth E. Jansen           enddo
136*59599516SKenneth E. Jansen
137*59599516SKenneth E. Jansen        enddo  ! quadrature point loop
138*59599516SKenneth E. Jansenc
139*59599516SKenneth E. Jansenc.... assemble the NABI vector
140*59599516SKenneth E. Jansenc
141*59599516SKenneth E. Jansen        call local (NABI,    rl,     ienb,   3,  'scatter ')
142*59599516SKenneth E. Jansenc
143*59599516SKenneth E. Jansenc     push the surf number which we have associated with boundary
144*59599516SKenneth E. JansenC     elements up to the global level in the array ndsurf
145*59599516SKenneth E. Jansenc
146*59599516SKenneth E. Jansen        do iel=1,npro
147*59599516SKenneth E. Jansen           if (iBCB(iel,2) .ne. 0) then
148*59599516SKenneth E. Jansen              iface = iBCB(iel,2)
149*59599516SKenneth E. Jansen              ndsurf(ienb(iel,1:nshlb))=iface
150*59599516SKenneth E. Jansen           endif
151*59599516SKenneth E. Jansen        enddo
152*59599516SKenneth E. Jansenc
153*59599516SKenneth E. Jansenc.... end
154*59599516SKenneth E. Jansenc
155*59599516SKenneth E. Jansen        return
156*59599516SKenneth E. Jansen        end
157*59599516SKenneth E. Jansen
158*59599516SKenneth E. Jansen
159*59599516SKenneth E. Jansen      subroutine AsBNASC ( x,       shpb,
160*59599516SKenneth E. Jansen     &                   ienb,  iBCB)
161*59599516SKenneth E. Jansenc
162*59599516SKenneth E. Jansenc----------------------------------------------------------------------
163*59599516SKenneth E. Jansenc
164*59599516SKenneth E. Jansenc This routine computes and assembles the data corresponding to the
165*59599516SKenneth E. Jansenc  boundary elements.
166*59599516SKenneth E. Jansenc
167*59599516SKenneth E. Jansenc Irene Vignon - copied from AsBNABI, Fall 2005.  (Fortran 90)
168*59599516SKenneth E. Jansenc----------------------------------------------------------------------
169*59599516SKenneth E. Jansenc
170*59599516SKenneth E. Jansen      use pvsQbi
171*59599516SKenneth E. Jansen        include "common.h"
172*59599516SKenneth E. Jansenc
173*59599516SKenneth E. Jansen        dimension xlb(npro,nenl,nsd),
174*59599516SKenneth E. Jansen     &            rl(npro,nshl),         WdetJb(npro),
175*59599516SKenneth E. Jansen     &            Wfactor(npro)
176*59599516SKenneth E. Jansen
177*59599516SKenneth E. Jansen        dimension x(numnp,nsd),
178*59599516SKenneth E. Jansen     &            shpb(nshl,ngaussb),    shglb(nsd,nshl,ngaussb),
179*59599516SKenneth E. Jansen     &            ienb(npro,nshl),
180*59599516SKenneth E. Jansen     &            iBCB(npro,ndiBCB)
181*59599516SKenneth E. Jansen
182*59599516SKenneth E. Jansen        dimension lnode(27),             sgn(npro,nshl),
183*59599516SKenneth E. Jansen     &            shpfun(npro,nshl),     shdrv(npro,nsd,nshl)
184*59599516SKenneth E. Jansen
185*59599516SKenneth E. Jansenc
186*59599516SKenneth E. Jansen        dimension dxdxib(npro,nsd,nsd),  temp(npro),
187*59599516SKenneth E. Jansen     &            temp1(npro),           temp2(npro),
188*59599516SKenneth E. Jansen     &            temp3(npro),
189*59599516SKenneth E. Jansen     &            v1(npro,nsd),          v2(npro,nsd)
190*59599516SKenneth E. Jansen
191*59599516SKenneth E. Jansenc
192*59599516SKenneth E. Jansenc.... get the matrix of mode signs for the hierarchic basis functions
193*59599516SKenneth E. Jansenc
194*59599516SKenneth E. Jansen        if (ipord .gt. 1) then
195*59599516SKenneth E. Jansen           call getsgn(ienb,sgn)
196*59599516SKenneth E. Jansen        endif
197*59599516SKenneth E. Jansenc
198*59599516SKenneth E. Jansenc.... gather the variables
199*59599516SKenneth E. Jansenc
200*59599516SKenneth E. Jansen        call localx(x,      xlb,    ienb,   nsd,    'gather  ')
201*59599516SKenneth E. Jansenc
202*59599516SKenneth E. Jansenc.... get the boundary element residuals
203*59599516SKenneth E. Jansenc
204*59599516SKenneth E. Jansen        rl  = zero
205*59599516SKenneth E. Jansenc
206*59599516SKenneth E. Jansenc.... compute the nodes which lie on the boundary (hierarchic)
207*59599516SKenneth E. Jansenc
208*59599516SKenneth E. Jansen        call getbnodes(lnode)
209*59599516SKenneth E. Jansenc
210*59599516SKenneth E. Jansenc.... compute the normal to the boundary. This is achieved by taking
211*59599516SKenneth E. Jansenc     the cross product of two vectors in the plane of the 2-d
212*59599516SKenneth E. Jansenc     boundary face.
213*59599516SKenneth E. Jansenc
214*59599516SKenneth E. Jansen  	if(lcsyst.eq.1) then	  ! set to curl into element all others out
215*59599516SKenneth E. Jansen  	   ipt2=2
216*59599516SKenneth E. Jansen  	   ipt3=3
217*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.2) then
218*59599516SKenneth E. Jansen  	   ipt2=4
219*59599516SKenneth E. Jansen  	   ipt3=2
220*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.3) then
221*59599516SKenneth E. Jansen  	   ipt2=3
222*59599516SKenneth E. Jansen  	   ipt3=2
223*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.4) then
224*59599516SKenneth E. Jansen  	   ipt2=2
225*59599516SKenneth E. Jansen  	   ipt3=4
226*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.5) then
227*59599516SKenneth E. Jansen  	   ipt2=4
228*59599516SKenneth E. Jansen  	   ipt3=2
229*59599516SKenneth E. Jansen  	elseif(lcsyst.eq.6) then
230*59599516SKenneth E. Jansen  	   ipt2=2
231*59599516SKenneth E. Jansen  	   ipt3=5
232*59599516SKenneth E. Jansen  	endif
233*59599516SKenneth E. Jansen  	v1 = xlb(:,ipt2,:) - xlb(:,1,:)
234*59599516SKenneth E. Jansen  	v2 = xlb(:,ipt3,:) - xlb(:,1,:)
235*59599516SKenneth E. Jansenc
236*59599516SKenneth E. Jansenc compute cross product
237*59599516SKenneth E. Jansenc
238*59599516SKenneth E. Jansen        temp1 = v1(:,2) * v2(:,3) - v2(:,2) * v1(:,3)
239*59599516SKenneth E. Jansen        temp2 = v2(:,1) * v1(:,3) - v1(:,1) * v2(:,3)
240*59599516SKenneth E. Jansen        temp3 = v1(:,1) * v2(:,2) - v2(:,1) * v1(:,2)
241*59599516SKenneth E. Jansenc
242*59599516SKenneth E. Jansenc mag is area for quads, twice area for tris
243*59599516SKenneth E. Jansenc
244*59599516SKenneth E. Jansen  	temp	   = one / sqrt ( temp1**2 + temp2**2 + temp3**2 )
245*59599516SKenneth E. Jansenc
246*59599516SKenneth E. Jansenc......here I only want the d Gamma, not the n_i
247*59599516SKenneth E. Jansen   	if (lcsyst .eq. 1) then
248*59599516SKenneth E. Jansen   	   Wfactor(:) = one / (four*temp(:))
249*59599516SKenneth E. Jansen   	elseif (lcsyst .eq. 2) then
250*59599516SKenneth E. Jansen   	   Wfactor(:) = one / (four*temp(:))
251*59599516SKenneth E. Jansen   	elseif (lcsyst .eq. 3) then
252*59599516SKenneth E. Jansen   	   Wfactor(:) = one / (two*temp(:))
253*59599516SKenneth E. Jansen   	elseif (lcsyst .eq. 4) then
254*59599516SKenneth E. Jansen   	   Wfactor(:) = one / (four*temp(:))
255*59599516SKenneth E. Jansen   	elseif (lcsyst .eq. 5) then
256*59599516SKenneth E. Jansen   	   Wfactor(:) = one / (four*temp(:))
257*59599516SKenneth E. Jansen   	elseif (lcsyst .eq. 6) then
258*59599516SKenneth E. Jansen   	   Wfactor(:) = one / (two*temp(:))
259*59599516SKenneth E. Jansen   	endif
260*59599516SKenneth E. Jansenc
261*59599516SKenneth E. Jansen        do iel=1,npro
262*59599516SKenneth E. Jansen           if (.not.btest(iBCB(iel,1),1)) then
263*59599516SKenneth E. Jansen              Wfactor(iel) = zero  ! we want zeros where we are not integrating
264*59599516SKenneth E. Jansen           endif
265*59599516SKenneth E. Jansen        enddo
266*59599516SKenneth E. Jansenc
267*59599516SKenneth E. Jansen        ngaussb = nintb(lcsyst)
268*59599516SKenneth E. Jansenc
269*59599516SKenneth E. Jansenc.... loop through the integration points
270*59599516SKenneth E. Jansenc
271*59599516SKenneth E. Jansen        do intp = 1, ngaussb
272*59599516SKenneth E. Jansen
273*59599516SKenneth E. Jansen           WdetJb(:) = Qwtb(lcsyst,intp) * Wfactor(:)
274*59599516SKenneth E. Jansenc
275*59599516SKenneth E. Jansenc.... get the hierarchic shape functions at this int point
276*59599516SKenneth E. Jansenc
277*59599516SKenneth E. Jansen           shglb=zero  ! protect debugger
278*59599516SKenneth E. Jansen           call getshpb(shpb,        shglb,        sgn,
279*59599516SKenneth E. Jansen     &              shpfun,       shdrv)
280*59599516SKenneth E. Jansen
281*59599516SKenneth E. Jansenc
282*59599516SKenneth E. Jansenc  Now lets calculate Integral N_(a:e) d Gamma
283*59599516SKenneth E. Jansenc
284*59599516SKenneth E. Jansen           do n = 1, nshlb
285*59599516SKenneth E. Jansen              nodlcl = lnode(n)
286*59599516SKenneth E. Jansen              rl(:,nodlcl) = rl(:,nodlcl) + shpfun(:,nodlcl) * WdetJb(:)
287*59599516SKenneth E. Jansen           enddo
288*59599516SKenneth E. Jansen
289*59599516SKenneth E. Jansen        enddo  ! quadrature point loop
290*59599516SKenneth E. Jansenc
291*59599516SKenneth E. Jansenc.... assemble the NASC vector
292*59599516SKenneth E. Jansenc
293*59599516SKenneth E. Jansen        call local (NASC,    rl,     ienb,   1,  'scatter ')
294*59599516SKenneth E. Jansenc
295*59599516SKenneth E. Jansenc     push the surf number which we have associated with boundary
296*59599516SKenneth E. JansenC     elements up to the global level in the array ndsurf
297*59599516SKenneth E. Jansenc
298*59599516SKenneth E. Jansen        do iel=1,npro
299*59599516SKenneth E. Jansen           if (iBCB(iel,2) .ne. 0) then
300*59599516SKenneth E. Jansen              iface = iBCB(iel,2)
301*59599516SKenneth E. Jansen              ndsurf(ienb(iel,1:nshlb))=iface
302*59599516SKenneth E. Jansen           endif
303*59599516SKenneth E. Jansen        enddo
304*59599516SKenneth E. Jansenc
305*59599516SKenneth E. Jansenc.... end
306*59599516SKenneth E. Jansenc
307*59599516SKenneth E. Jansen        return
308*59599516SKenneth E. Jansen        end
309