CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
slavintpoints.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine slavintpoints (ntie, itietri, ipkon, kon, lakon, straight, nintpoint, koncont, co, vold, xo, yo, zo, x, y, z, nx, ny, nz, islavsurf, islavnode, nslavnode, imastop, mi, ncont, ipe, ime, pslavsurf, i, l, ntri)
 

Function/Subroutine Documentation

◆ slavintpoints()

subroutine slavintpoints ( integer  ntie,
integer, dimension(2,ntie)  itietri,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
real*8, dimension(16,*)  straight,
integer  nintpoint,
integer, dimension(4,*)  koncont,
real*8, dimension(3,*)  co,
real*8, dimension(0:mi(2),*)  vold,
real*8, dimension(*)  xo,
real*8, dimension(*)  yo,
real*8, dimension(*)  zo,
real*8, dimension(*)  x,
real*8, dimension(*)  y,
real*8, dimension(*)  z,
integer, dimension(*)  nx,
integer, dimension(*)  ny,
integer, dimension(*)  nz,
integer, dimension(2,*)  islavsurf,
integer, dimension(*)  islavnode,
integer, dimension(ntie+1)  nslavnode,
integer, dimension(3,*)  imastop,
integer, dimension(*)  mi,
integer  ncont,
integer, dimension(*)  ipe,
integer, dimension(4,*)  ime,
real*8, dimension(3,*)  pslavsurf,
integer  i,
integer  l,
integer  ntri 
)
30 !
31 ! Author: Li, Yang; Rakotonanahary, Samoela; Sitzmann,Saskia
32 !
33  implicit none
34 !
35  character*8 lakon(*)
36 !
37  integer ntie,nintpoint,imastop(3,*),ncont,itietri(2,ntie),
38  & ipkon(*),kon(*),koncont(4,*),node,neigh(10),iflag,kneigh,i,
39  & j,k,l,ii,itri,nx(*),ny(*),nz(*),ifreeintersec,nelemm,jfacem,
40  & indexe,nnodelem,nope,islavsurf(2,*),islavnode(*),
41  & nslavnode(ntie+1),ifaces,nelems,jfaces,mi(*),
42  & m,nopes,konl(20),id,imface(8),ntria,
43  & imfacecorner(8,8),line,iactiveline(3,3*ncont),
44  & icoveredmelem(3*ncont),nactiveline,ipe(*),ime(4,*),k1,j1,
45  & ncoveredmelem,ntri,nintpfirst,nodem(8),nodem2(8),
46  & il,nodel,getnodel,ifacem,idummy,nopemm,nmp,k2,j2
47 !
48  real*8 straight(16,*),co(3,*),vold(0:mi(2),*),
49  & xo(*),yo(*),zo(*),x(*),y(*),z(*),
50  & xl2m(3,8),xl2s(3,8),xl2m2(3,8),
51  & pmiddle(3),xl2sr(3,8),xl2sp(3,8),xl2s2(3,8),
52  & dd,xns(3,8),areaslav,al,xn(3),slavstraight(36),
53  & err2,dist,distmin,pslavsurf(3,*),err,xquad(2,8),
54  & xtri(2,6),xi,et,xsj2(3),xs2(3,2),shp2(7,8),anglesm
55 !
56  data iflag /2/
57 !
58  data xquad /-1,-1,
59  & 1,-1,
60  & 1,1,
61  & -1,1,
62  & 0,-1,
63  & 1,0,
64  & 0,1,
65  & -1,0/
66 !
67  data xtri /0,0,
68  & 1,0,
69  & 0,1,
70  & 0.5,0,
71  & 0.5,0.5,
72  & 0,0.5/
73 !
74  kneigh=1
75  err=0.1
76  areaslav=0.0
77  nintpfirst=nintpoint
78  islavsurf(2,l)=nintpoint
79 !
80 ! Research of the contact integration points
81 !
82  ifaces=islavsurf(1,l)
83  nelems=int(ifaces/10)
84  jfaces=ifaces-nelems*10
85 !
86 ! get nope,nopes
87 !
88  call faceinfo(nelems,jfaces,lakon,nope,nopes,idummy)
89 !
90 ! actual position of the nodes belonging to the
91 ! slave surface
92 !
93  do j=1,nope
94  konl(j)=kon(ipkon(nelems)+j)
95  enddo
96 !
97  do m=1,nopes
98  do j=1,3
99  nodel=getnodel(m,jfaces,nope)
100  xl2s(j,m)=co(j,konl(nodel))+
101  & vold(j,konl(nodel))
102  enddo
103  enddo
104 !
105 ! slightly reducing the size of the slave surface in
106 ! an aleatoric way
107 !
108  do j=1,3
109  pmiddle(j)=0.d0
110  do m=1,nopes
111  pmiddle(j)=pmiddle(j)+xl2s(j,m)
112  enddo
113  pmiddle(j)=pmiddle(j)/nopes
114  enddo
115  do j=1,3
116  do m=1,nopes
117  xl2sr(j,m)=xl2s(j,m)-0.5*err*(xl2s(j,m)-pmiddle(j))
118  enddo
119  enddo
120 !
121 ! Resort vertices for quadratic elements
122 !
123  if((nopes.eq.3).or.(nopes.eq.4))then
124  do j=1,nopes
125  do k=1,3
126  xl2s2(k,j)=xl2s(k,j)
127  enddo
128  enddo
129  else
130  do j=1,int(nopes/2.0)
131  do k=1,3
132  xl2s2(k,2*j-1)=xl2s(k,j)
133  xl2s2(k,2*j)=xl2s(k,(int(nopes/2.0))+j)
134  enddo
135  enddo
136  endif
137 !
138 ! calculate the mean normal vector on the Slave Surface
139 ! +
140 ! determine the equations of the triangle/quadrilateral
141 ! (mean)plane and of the planes boardering the
142 ! triangle/quadrilateral
143 !
144  if(nopes.eq.3) then
145  call straighteq3d(xl2s2,slavstraight)
146  do k=1,3
147  xn(k)=slavstraight(4*nopes+k)
148  enddo
149  else
150  do k=1,3
151  xn(k)=0.d0
152  enddo
153 !
154  do m=1,nopes
155  if((nopes.eq.4).or.(nopes.eq.8))then
156  xi=xquad(1,m)
157  et=xquad(2,m)
158  else
159  xi=xtri(1,m)
160  et=xtri(2,m)
161  endif
162  if(nopes.eq.8)then
163  call shape8q(xi,et,xl2s,xsj2,xs2,shp2,iflag)
164  elseif(nopes.eq.4)then
165  call shape4q(xi,et,xl2s,xsj2,xs2,shp2,iflag)
166  elseif(nopes.eq.6)then
167  call shape6tri(xi,et,xl2s,xsj2,xs2,shp2,iflag)
168  else
169  call shape3tri(xi,et,xl2s,xsj2,xs2,shp2,iflag)
170  endif
171  dd=dsqrt(xsj2(1)*xsj2(1)+xsj2(2)*xsj2(2)
172  & +xsj2(3)*xsj2(3))
173  xsj2(1)=xsj2(1)/dd
174  xsj2(2)=xsj2(2)/dd
175  xsj2(3)=xsj2(3)/dd
176 !
177  do k=1,3
178  xn(k)=xn(k)
179  & +xsj2(k)
180  enddo
181  enddo
182 !
183 ! normalizing the mean normal on the Slave surface
184 !
185  dd=dsqrt(xn(1)**2+xn(2)**2+xn(3)**2)
186  do k=1,3
187  xn(k)=xn(k)/dd
188  enddo
189  call approxplane(xl2s2,slavstraight,xn,nopes)
190  endif
191 !
192 ! Project slave nodes to meanplane, needed for Sutherland-Hodgman
193 !
194  do j=1,nopes
195  al=-xn(1)*xl2s2(1,j)-xn(2)*
196  & xl2s2(2,j)-xn(3)*xl2s2(3,j)-
197  & slavstraight(nopes*4+4)
198  if(nopes.ne.3)then
199  do k=1,3
200  xl2sp(k,j)=xl2s2(k,j)+al*xn(k)
201  enddo
202  else
203  do k=1,3
204  xl2sp(k,j)=xl2s2(k,j)
205  enddo
206  endif
207  enddo
208 !
209 ! determine the triangles corresponding to the corner
210 ! nodes
211 !
212  ntria=0
213  do j=1,8
214  imface(j)=0
215  do k=1,8
216  imfacecorner(j,k)=0
217  enddo
218  enddo
219 !
220  do j=1,3
221  do m=1,nopes
222  xl2sr(j,m)=xl2s(j,m)-2*err*(xl2s(j,m)-pmiddle(j))
223  enddo
224  enddo
225  distmin=1.1
226  do j=1,nopes
227  call neartriangle(xl2sr(1,j),xn,xo,yo,zo,x,y,z,nx,ny,nz,
228  & ntri,neigh,kneigh,itietri,ntie,straight,imastop,itri,i)
229  nodel=getnodel(j,jfaces,nope)
230  node=konl(nodel)
231  if(itri.eq.0) then
232  cycle
233  endif
234  dist=-(straight(13,itri)*xl2sr(1,j)+
235  & straight(14,itri)*xl2sr(2,j)+
236  & straight(15,itri)*xl2sr(3,j)+
237  & straight(16,itri))/
238  & (straight(13,itri)*xn(1)+
239  & straight(14,itri)*xn(2)+
240  & straight(15,itri)*xn(3))
241  if(dist.lt.distmin)distmin=dist
242  ifacem=koncont(4,itri)
243 !
244  call nident(imface,ifacem,ntria,id)
245  if(id.gt.0) then
246  if(imface(id).eq.ifacem) then
247  imfacecorner(j,id)=1
248  cycle
249  endif
250  endif
251 !
252 ! triangle was not covered yet: add to stack
253 !
254 ! angle criteria
255 !
256  anglesm=xn(1)*straight(13,itri)
257  & +xn(2)*straight(14,itri)
258  & +xn(3)*straight(15,itri)
259 !
260  if(anglesm.lt.-0.7)then
261 ! angle between surface normals between 135 and 225 degree
262 ! angle between surfaces between 0 and 45 degree
263  ntria=ntria+1
264  do k=ntria,id+2,-1
265  imface(k)=imface(k-1)
266  do m=1,j-1
267  imfacecorner(m,k)=imfacecorner(m,k-1)
268  enddo
269  enddo
270  imface(id+1)=ifacem
271  imfacecorner(j,id+1)=1
272  do m=1,j-1
273  imfacecorner(m,id+1)=0
274  enddo
275  endif
276  enddo
277 !
278  nactiveline=0
279  ifreeintersec=0
280 !
281 ! treating the corner elements first
282 !
283  ncoveredmelem=0
284  do j=1,ntria
285  ifacem=imface(j)
286  nelemm=int(ifacem/10.d0)
287  jfacem=ifacem-10*nelemm
288 !
289 ! add master element to covered stack
290 !
291  call nident(icoveredmelem,nelemm,ncoveredmelem,id)
292  if(id.ne.0 .and. icoveredmelem(id).eq.nelemm)then
293 ! master element was already threated
294  cycle
295  else
296 ! add master element to covered elements
297  ncoveredmelem=ncoveredmelem+1
298  do ii=ncoveredmelem,id+2,-1
299  icoveredmelem(ii)=icoveredmelem(ii-1)
300  enddo
301  icoveredmelem(id+1)=nelemm
302  endif
303  call faceinfo(nelemm,jfacem,lakon,nopemm,
304  & nnodelem,idummy)
305 !
306 ! determining the nodes of the face
307 !
308  do j1=1,nopemm
309  konl(j1)=kon(ipkon(nelemm)+j1)
310  enddo
311  do k1=1,nnodelem
312  nodel=getnodel(k1,jfacem,nopemm)
313  nodem(k1)=konl(nodel)
314  do j1=1,3
315  xl2m(j1,k1)=co(j1,konl(nodel))+
316  & vold(j1,konl(nodel))
317  enddo
318  enddo
319  dd=dsqrt(xn(1)**2+xn(2)**2+xn(3)**2)
320 !
321 ! divide master element into konvex subelements
322 !
323  if((nnodelem.eq.3).or.(nnodelem.eq.4))then
324 !
325 ! no loop needed
326 !
327  nmp=nnodelem
328  do k2=1,nnodelem
329  nodem2(k2)=nodem(k2)
330  do j2=1,3
331  xl2m2(j2,k2)=xl2m(j2,k2)
332  enddo
333  enddo
334  call treatmasterface(
335  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
336  & ipe,ime,iactiveline,nactiveline,
337  & ifreeintersec,ifacem,
338  & nintpoint,pslavsurf,
339  & xl2m,nnodelem,xl2m2,nmp,nodem2,
340  & areaslav)
341 !
342  elseif(nnodelem.eq.6)then
343 !
344 ! tri6 surface is divided into 4 tri3 surfaces
345 !
346 ! 1. triangle
347 !
348  nmp=3
349  nodem2(1)=nodem(1)
350  nodem2(2)=nodem(4)
351  nodem2(3)=nodem(6)
352  do j2=1,3
353  xl2m2(j2,1)=xl2m(j2,1)
354  enddo
355  do j2=1,3
356  xl2m2(j2,2)=xl2m(j2,4)
357  enddo
358  do j2=1,3
359  xl2m2(j2,3)=xl2m(j2,6)
360  enddo
361  call treatmasterface(
362  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
363  & ipe,ime,iactiveline,nactiveline,
364  & ifreeintersec,ifacem,
365  & nintpoint,pslavsurf,
366  & xl2m,nnodelem,xl2m2,nmp,nodem2,
367  & areaslav)
368 !
369 ! 2. triangle
370 !
371  nmp=3
372  nodem2(1)=nodem(4)
373  nodem2(2)=nodem(2)
374  nodem2(3)=nodem(5)
375  do j2=1,3
376  xl2m2(j2,1)=xl2m(j2,4)
377  enddo
378  do j2=1,3
379  xl2m2(j2,2)=xl2m(j2,2)
380  enddo
381  do j2=1,3
382  xl2m2(j2,3)=xl2m(j2,5)
383  enddo
384  call treatmasterface(
385  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
386  & ipe,ime,iactiveline,nactiveline,
387  & ifreeintersec,ifacem,
388  & nintpoint,pslavsurf,
389  & xl2m,nnodelem,xl2m2,nmp,nodem2,
390  & areaslav)
391 !
392 ! 3. triangle
393 !
394  nmp=3
395  nodem2(1)=nodem(5)
396  nodem2(2)=nodem(3)
397  nodem2(3)=nodem(6)
398  do j2=1,3
399  xl2m2(j2,1)=xl2m(j2,5)
400  enddo
401  do j2=1,3
402  xl2m2(j2,2)=xl2m(j2,3)
403  enddo
404  do j2=1,3
405  xl2m2(j2,3)=xl2m(j2,6)
406  enddo
407  call treatmasterface(
408  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
409  & ipe,ime,iactiveline,nactiveline,
410  & ifreeintersec,ifacem,
411  & nintpoint,pslavsurf,
412  & xl2m,nnodelem,xl2m2,nmp,nodem2,
413  & areaslav)
414 !
415 ! 4. triangle
416 !
417  nmp=3
418  nodem2(1)=nodem(4)
419  nodem2(2)=nodem(5)
420  nodem2(3)=nodem(6)
421  do j2=1,3
422  xl2m2(j2,1)=xl2m(j2,4)
423  enddo
424  do j2=1,3
425  xl2m2(j2,2)=xl2m(j2,5)
426  enddo
427  do j2=1,3
428  xl2m2(j2,3)=xl2m(j2,6)
429  enddo
430  call treatmasterface(
431  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
432  & ipe,ime,iactiveline,nactiveline,
433  & ifreeintersec,ifacem,
434  & nintpoint,pslavsurf,
435  & xl2m,nnodelem,xl2m2,nmp,nodem2,
436  & areaslav)
437 !
438  elseif(nnodelem.eq.8)then
439 !
440 ! quad8 surface is divided into 4 tri3 + 1 quad4 surfaces
441 !
442 ! 1. triangle
443 !
444  nmp=3
445  nodem2(1)=nodem(1)
446  nodem2(2)=nodem(5)
447  nodem2(3)=nodem(8)
448  do j2=1,3
449  xl2m2(j2,1)=xl2m(j2,1)
450  enddo
451  do j2=1,3
452  xl2m2(j2,2)=xl2m(j2,5)
453  enddo
454  do j2=1,3
455  xl2m2(j2,3)=xl2m(j2,8)
456  enddo
457  call treatmasterface(
458  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
459  & ipe,ime,iactiveline,nactiveline,
460  & ifreeintersec,ifacem,
461  & nintpoint,pslavsurf,
462  & xl2m,nnodelem,xl2m2,nmp,nodem2,
463  & areaslav)
464 !
465 ! 2. triangle
466 !
467  nmp=3
468  nodem2(1)=nodem(5)
469  nodem2(2)=nodem(2)
470  nodem2(3)=nodem(6)
471  do j2=1,3
472  xl2m2(j2,1)=xl2m(j2,5)
473  enddo
474  do j2=1,3
475  xl2m2(j2,2)=xl2m(j2,2)
476  enddo
477  do j2=1,3
478  xl2m2(j2,3)=xl2m(j2,6)
479  enddo
480  call treatmasterface(
481  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
482  & ipe,ime,iactiveline,nactiveline,
483  & ifreeintersec,ifacem,
484  & nintpoint,pslavsurf,
485  & xl2m,nnodelem,xl2m2,nmp,nodem2,
486  & areaslav)
487 !
488 ! 3. triangle
489 !
490  nmp=3
491  nodem2(1)=nodem(6)
492  nodem2(2)=nodem(3)
493  nodem2(3)=nodem(7)
494  do j2=1,3
495  xl2m2(j2,1)=xl2m(j2,6)
496  enddo
497  do j2=1,3
498  xl2m2(j2,2)=xl2m(j2,3)
499  enddo
500  do j2=1,3
501  xl2m2(j2,3)=xl2m(j2,7)
502  enddo
503  call treatmasterface(
504  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
505  & ipe,ime,iactiveline,nactiveline,
506  & ifreeintersec,ifacem,
507  & nintpoint,pslavsurf,
508  & xl2m,nnodelem,xl2m2,nmp,nodem2,
509  & areaslav)
510 !
511 ! 4. triangle
512 !
513  nmp=3
514  nodem2(1)=nodem(7)
515  nodem2(2)=nodem(4)
516  nodem2(3)=nodem(8)
517  do j2=1,3
518  xl2m2(j2,1)=xl2m(j2,7)
519  enddo
520  do j2=1,3
521  xl2m2(j2,2)=xl2m(j2,4)
522  enddo
523  do j2=1,3
524  xl2m2(j2,3)=xl2m(j2,8)
525  enddo
526  call treatmasterface(
527  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
528  & ipe,ime,iactiveline,nactiveline,
529  & ifreeintersec,ifacem,
530  & nintpoint,pslavsurf,
531  & xl2m,nnodelem,xl2m2,nmp,nodem2,
532  & areaslav)
533 !
534 ! quad
535 !
536  nmp=4
537  nodem2(1)=nodem(5)
538  nodem2(2)=nodem(6)
539  nodem2(3)=nodem(7)
540  nodem2(4)=nodem(8)
541  do j2=1,3
542  xl2m2(j2,1)=xl2m(j2,5)
543  enddo
544  do j2=1,3
545  xl2m2(j2,2)=xl2m(j2,6)
546  enddo
547  do j2=1,3
548  xl2m2(j2,3)=xl2m(j2,7)
549  enddo
550  do j2=1,3
551  xl2m2(j2,4)=xl2m(j2,8)
552  enddo
553  call treatmasterface(
554  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
555  & ipe,ime,iactiveline,nactiveline,
556  & ifreeintersec,ifacem,
557  & nintpoint,pslavsurf,
558  & xl2m,nnodelem,xl2m2,nmp,nodem2,
559  & areaslav)
560  endif
561  enddo
562 !
563 ! corners of the Slave surface have already been treated
564 !
565  do j=1,nopes
566  imfacecorner(j,1)=0
567  enddo
568 !
569 ! retrieving all triangles by neighborhood search
570 !
571  do
572  line=iactiveline(1,1)
573  if(nactiveline.eq.0) exit
574  if(koncont(4,ime(2,line)).eq.iactiveline(2,1)) then
575  itri=imastop(ime(3,line),ime(2,line))
576  else
577  itri=ime(2,line)
578  endif
579 !
580 ! check whether still in contact tie
581 !
582  if((itri.gt.itietri(2,i)).or.(itri.lt.itietri(1,i)))then
583  if(itri.ne.0)then
584  endif
585  itri=0
586  endif
587 
588  if(itri.eq.0) then
589  nactiveline=nactiveline-1
590  do il=1,nactiveline
591  do k=1,3
592  iactiveline(k,il)=iactiveline(k,il+1)
593  enddo
594  enddo
595  cycle
596  endif
597 
598 !
599  ifacem=koncont(4,itri)
600  nelemm=int(koncont(4,itri)/10.d0)
601  jfacem=koncont(4,itri)-10*nelemm
602 !
603 ! add master element to covered stack
604 !
605  call nident(icoveredmelem,nelemm,ncoveredmelem,id)
606  if(id .gt. 0 .and. icoveredmelem(id).eq.nelemm)then
607 !
608 ! master elment was already threated
609 !
610  nactiveline=nactiveline-1
611  do il=1,nactiveline
612  do k=1,3
613  iactiveline(k,il)=iactiveline(k,il+1)
614  enddo
615  enddo
616  cycle
617  else
618 !
619 ! add master element to covered elements
620 !
621  ncoveredmelem=ncoveredmelem+1
622  do ii=ncoveredmelem,id+2,-1
623  icoveredmelem(ii)=icoveredmelem(ii-1)
624  enddo
625  icoveredmelem(id+1)=nelemm
626  endif
627  indexe=ipkon(nelemm)
628  call faceinfo(nelemm,jfacem,lakon,nopemm,
629  & nnodelem,idummy)
630 !
631 ! determining the nodes of the face
632 !
633  do j1=1,nopemm
634  konl(j1)=kon(ipkon(nelemm)+j1)
635  enddo
636  do k1=1,nnodelem
637  nodel=getnodel(k1,jfacem,nopemm)
638  nodem(k1)=konl(nodel)
639  do j1=1,3
640  xl2m(j1,k1)=co(j1,konl(nodel))+
641  & vold(j1,konl(nodel))
642  enddo
643  enddo
644 !
645 ! divide master element into konvex subelements
646 !
647  if((nnodelem.eq.3).or.(nnodelem.eq.4))then
648 !
649 ! no loop needed
650 !
651  nmp=nnodelem
652  do k2=1,nnodelem
653  nodem2(k2)=nodem(k2)
654  do j2=1,3
655  xl2m2(j2,k2)=xl2m(j2,k2)
656  enddo
657  enddo
658  call treatmasterface(
659  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
660  & ipe,ime,iactiveline,nactiveline,
661  & ifreeintersec,ifacem,
662  & nintpoint,pslavsurf,
663  & xl2m,nnodelem,xl2m2,nmp,nodem2,
664  & areaslav)
665 !
666  elseif(nnodelem.eq.6)then
667 !
668 ! tri6 surface is divided into 4 tri3 surfaces
669 !
670 ! 1. triangle
671 !
672  nmp=3
673  nodem2(1)=nodem(1)
674  nodem2(2)=nodem(4)
675  nodem2(3)=nodem(6)
676  do j2=1,3
677  xl2m2(j2,1)=xl2m(j2,1)
678  enddo
679  do j2=1,3
680  xl2m2(j2,2)=xl2m(j2,4)
681  enddo
682  do j2=1,3
683  xl2m2(j2,3)=xl2m(j2,6)
684  enddo
685  call treatmasterface(
686  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
687  & ipe,ime,iactiveline,nactiveline,
688  & ifreeintersec,ifacem,
689  & nintpoint,pslavsurf,
690  & xl2m,nnodelem,xl2m2,nmp,nodem2,
691  & areaslav)
692 !
693 ! 2. triangle
694 !
695  nmp=3
696  nodem2(1)=nodem(4)
697  nodem2(2)=nodem(2)
698  nodem2(3)=nodem(5)
699  do j2=1,3
700  xl2m2(j2,1)=xl2m(j2,4)
701  enddo
702  do j2=1,3
703  xl2m2(j2,2)=xl2m(j2,2)
704  enddo
705  do j2=1,3
706  xl2m2(j2,3)=xl2m(j2,5)
707  enddo
708  call treatmasterface(
709  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
710  & ipe,ime,iactiveline,nactiveline,
711  & ifreeintersec,ifacem,
712  & nintpoint,pslavsurf,
713  & xl2m,nnodelem,xl2m2,nmp,nodem2,
714  & areaslav)
715 !
716 ! 3. triangle
717 !
718  nmp=3
719  nodem2(1)=nodem(5)
720  nodem2(2)=nodem(3)
721  nodem2(3)=nodem(6)
722  do j2=1,3
723  xl2m2(j2,1)=xl2m(j2,5)
724  enddo
725  do j2=1,3
726  xl2m2(j2,2)=xl2m(j2,3)
727  enddo
728  do j2=1,3
729  xl2m2(j2,3)=xl2m(j2,6)
730  enddo
731  call treatmasterface(
732  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
733  & ipe,ime,iactiveline,nactiveline,
734  & ifreeintersec,ifacem,
735  & nintpoint,pslavsurf,
736  & xl2m,nnodelem,xl2m2,nmp,nodem2,
737  & areaslav)
738 !
739 ! 4. triangle
740 !
741  nmp=3
742  nodem2(1)=nodem(4)
743  nodem2(2)=nodem(5)
744  nodem2(3)=nodem(6)
745  do j2=1,3
746  xl2m2(j2,1)=xl2m(j2,4)
747  enddo
748  do j2=1,3
749  xl2m2(j2,2)=xl2m(j2,5)
750  enddo
751  do j2=1,3
752  xl2m2(j2,3)=xl2m(j2,6)
753  enddo
754  call treatmasterface(
755  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
756  & ipe,ime,iactiveline,nactiveline,
757  & ifreeintersec,ifacem,
758  & nintpoint,pslavsurf,
759  & xl2m,nnodelem,xl2m2,nmp,nodem2,
760  & areaslav)
761 !
762  elseif(nnodelem.eq.8)then
763 !
764 ! quad8 surface is divided into 4 tri3 + 1 quad4 surfaces
765 !
766 ! 1. triangle
767 !
768  nmp=3
769  nodem2(1)=nodem(1)
770  nodem2(2)=nodem(5)
771  nodem2(3)=nodem(8)
772  do j2=1,3
773  xl2m2(j2,1)=xl2m(j2,1)
774  enddo
775  do j2=1,3
776  xl2m2(j2,2)=xl2m(j2,5)
777  enddo
778  do j2=1,3
779  xl2m2(j2,3)=xl2m(j2,8)
780  enddo
781  call treatmasterface(
782  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
783  & ipe,ime,iactiveline,nactiveline,
784  & ifreeintersec,ifacem,
785  & nintpoint,pslavsurf,
786  & xl2m,nnodelem,xl2m2,nmp,nodem2,
787  & areaslav)
788 !
789 ! 2. triangle
790 !
791  nmp=3
792  nodem2(1)=nodem(5)
793  nodem2(2)=nodem(2)
794  nodem2(3)=nodem(6)
795  do j2=1,3
796  xl2m2(j2,1)=xl2m(j2,5)
797  enddo
798  do j2=1,3
799  xl2m2(j2,2)=xl2m(j2,2)
800  enddo
801  do j2=1,3
802  xl2m2(j2,3)=xl2m(j2,6)
803  enddo
804  call treatmasterface(
805  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
806  & ipe,ime,iactiveline,nactiveline,
807  & ifreeintersec,ifacem,
808  & nintpoint,pslavsurf,
809  & xl2m,nnodelem,xl2m2,nmp,nodem2,
810  & areaslav)
811 !
812 ! 3. triangle
813 !
814  nmp=3
815  nodem2(1)=nodem(6)
816  nodem2(2)=nodem(3)
817  nodem2(3)=nodem(7)
818  do j2=1,3
819  xl2m2(j2,1)=xl2m(j2,6)
820  enddo
821  do j2=1,3
822  xl2m2(j2,2)=xl2m(j2,3)
823  enddo
824  do j2=1,3
825  xl2m2(j2,3)=xl2m(j2,7)
826  enddo
827  call treatmasterface(
828  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
829  & ipe,ime,iactiveline,nactiveline,
830  & ifreeintersec,ifacem,
831  & nintpoint,pslavsurf,
832  & xl2m,nnodelem,xl2m2,nmp,nodem2,
833  & areaslav)
834 !
835 ! 4. triangle
836 !
837  nmp=3
838  nodem2(1)=nodem(7)
839  nodem2(2)=nodem(4)
840  nodem2(3)=nodem(8)
841  do j2=1,3
842  xl2m2(j2,1)=xl2m(j2,7)
843  enddo
844  do j2=1,3
845  xl2m2(j2,2)=xl2m(j2,4)
846  enddo
847  do j2=1,3
848  xl2m2(j2,3)=xl2m(j2,8)
849  enddo
850  call treatmasterface(
851  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
852  & ipe,ime,iactiveline,nactiveline,
853  & ifreeintersec,ifacem,
854  & nintpoint,pslavsurf,
855  & xl2m,nnodelem,xl2m2,nmp,nodem2,
856  & areaslav)
857 !
858 ! quad
859 !
860  nmp=4
861  nodem2(1)=nodem(5)
862  nodem2(2)=nodem(6)
863  nodem2(3)=nodem(7)
864  nodem2(4)=nodem(8)
865  do j2=1,3
866  xl2m2(j2,1)=xl2m(j2,5)
867  enddo
868  do j2=1,3
869  xl2m2(j2,2)=xl2m(j2,6)
870  enddo
871  do j2=1,3
872  xl2m2(j2,3)=xl2m(j2,7)
873  enddo
874  do j2=1,3
875  xl2m2(j2,4)=xl2m(j2,8)
876  enddo
877  call treatmasterface(
878  & nopes,slavstraight,xn,xns,xl2s,xl2sp,
879  & ipe,ime,iactiveline,nactiveline,
880  & ifreeintersec,ifacem,
881  & nintpoint,pslavsurf,
882  & xl2m,nnodelem,xl2m2,nmp,nodem2,
883  & areaslav)
884  endif
885  enddo
886  islavsurf(2,l+1)=nintpoint
887 !
888  return
subroutine shape8q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape8q.f:20
subroutine shape3tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape3tri.f:20
subroutine faceinfo(nelem, jface, lakon, nope, nopes, mint2d)
Definition: faceinfo.f:24
integer function getnodel(ii, jface, nope)
Definition: getnodel.f:24
static double * dist
Definition: radflowload.c:42
subroutine shape4q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape4q.f:20
subroutine neartriangle(p, xn, xo, yo, zo, x, y, z, nx, ny, nz, n, neigh, kneigh, itietri, ntie, straight, imastop, itri, itie)
Definition: neartriangle.f:22
subroutine treatmasterface(nopes, slavstraight, xn, xns, xl2s, xl2sp, ipe, ime, iactiveline, nactiveline, ifreeintersec, nelemm, nintpoint, pslavsurf, xl2m, nnodelem, xl2m2, nmp, nodem, areaslav)
Definition: treatmasterface.f:27
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine straighteq3d(col, straight)
Definition: straighteq3d.f:20
subroutine shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
subroutine approxplane(col, straight, xn, nopes)
Definition: approxplane.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)