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

Go to the source code of this file.

Functions/Subroutines

subroutine treatmasterface (nopes, slavstraight, xn, xns, xl2s, xl2sp, ipe, ime, iactiveline, nactiveline, ifreeintersec, nelemm, nintpoint, pslavsurf, xl2m, nnodelem, xl2m2, nmp, nodem, areaslav)
 

Function/Subroutine Documentation

◆ treatmasterface()

subroutine treatmasterface ( integer, intent(in)  nopes,
real*8, dimension(36), intent(in)  slavstraight,
real*8, dimension(3), intent(in)  xn,
real*8, dimension(3,8), intent(in)  xns,
real*8, dimension(3,*), intent(in)  xl2s,
real*8, dimension(3,*), intent(in)  xl2sp,
integer, dimension(*), intent(in)  ipe,
integer, dimension(4,*), intent(in)  ime,
integer, dimension(3,*), intent(in)  iactiveline,
integer, intent(inout)  nactiveline,
integer, intent(in)  ifreeintersec,
integer, intent(in)  nelemm,
integer, intent(inout)  nintpoint,
real*8, dimension(3,*), intent(inout)  pslavsurf,
real*8, dimension(3,8), intent(in)  xl2m,
integer, intent(in)  nnodelem,
real*8, dimension(3,8), intent(in)  xl2m2,
integer, intent(in)  nmp,
integer, dimension(*), intent(in)  nodem,
real*8, intent(inout)  areaslav 
)
27 !
28 ! Author: Saskia Sitzmann
29 !
30  implicit none
31 !
32  integer nvertex,nopes,ipe(*),ime(4,*),iactiveline(3,*),
33  & nactiveline,ifreeintersec,nmp,i,j,k,nintpoint,
34  & nnodelem,nodem(*),modf,nelemm,k_max,nipold
35 !
36  real*8 pvertex(3,13),slavstraight(36),xn(3),xilm,etlm,xnl(3),
37  & xl2s(3,*),p1(2),p2(2),pslavsurf(3,*),xil,etl,
38  & area,xl2m(3,8),xl2m2(3,8),al,err,xns(3,8),
39  & xl2sp(3,*),xl2mp(3,8),cgp(3),pm(3),ps(3),xit(3),etat(3),areaslav
40 !
41  intent(in) nopes,slavstraight,xn,xns,xl2s,xl2sp,
42  & ipe,ime,iactiveline,ifreeintersec,nelemm,
43  & xl2m,nnodelem,xl2m2,nmp,nodem
44 !
45  intent(inout) pslavsurf,nintpoint,areaslav,nactiveline
46 !
47  include "gauss.f"
48 !
49  err=1.d-6
50  nvertex=0
51  nipold=nintpoint
52 !
53 ! Project master nodes to meanplane, needed for Sutherland-Hodgman
54 !
55  do j=1, nmp
56  al=-xn(1)*xl2m2(1,j)-xn(2)*
57  & xl2m2(2,j)-xn(3)*
58  & xl2m2(3,j)-slavstraight(nopes*4+4)
59  do k=1,3
60  xl2mp(k,j)= xl2m2(k,j)+al*xn(k)
61  enddo
62  enddo
63 !
64 ! call Sutherland-Hodgman Algo
65 !
66  call sutherland_hodgman(nopes,xn,xl2sp,xl2mp,nodem,
67  & ipe,ime,iactiveline,nactiveline,
68  & ifreeintersec,nelemm,nmp,
69  & nvertex,pvertex)
70 !
71 !
72  do k=1,3
73  cgp(k)=0.0
74  enddo
75  if(nvertex.lt.3) return
76 !
77  if(nvertex==3)then
78  do k=1,3
79  cgp(k)=pvertex(k,nvertex)
80  enddo
81  nvertex=nvertex-1
82  k_max=1
83  else
84  do i=1,nvertex
85  do k=1,3
86  cgp(k)=cgp(k)+pvertex(k,i)/nvertex
87  enddo
88  enddo
89  k_max=nvertex
90  endif
91 !
92 ! Project center point back on slave face
93 !
94  call attachline(xl2s,cgp,nopes,xit(3),etat(3),xn)
95 !
96 ! generating integration points on the slave surface S
97 !
98  do k=1,k_max
99 !
100 ! Project back on slave surface
101 !
102  call attachline(xl2s,pvertex(1:3,modf(nvertex,k)),
103  & nopes,xit(1),etat(1),xn)
104  call attachline(xl2s,pvertex(1:3,modf(nvertex,k+1)),
105  & nopes,xit(2),etat(2),xn)
106 !
107  p1(1)=xit(1)-xit(3)
108  p1(2)=etat(1)-etat(3)
109 !
110  p2(1)=xit(2)-xit(3)
111  p2(2)=etat(2)-etat(3)
112 !
113  area=dabs(p1(1)*p2(2)-p2(1)*p1(2))/2.d0
114 !
115  if(area.lt.1.e-4) cycle
116  if(nopes.eq.4.and.areaslav+area-4.0.gt.1.e-3
117  & .and.nactiveline.gt.0)then
118  nactiveline=0
119  return
120  endif
121  if(nopes.eq.3.and.areaslav+area-0.5.gt.1.e-4
122  & .and.nactiveline.gt.0)then
123  nactiveline=0
124  return
125  endif
126  areaslav=areaslav+area
127 !
128 ! 7 points scheme
129 !
130  do i=1,7
131  xil= xit(3)*gauss2d6(1,i)+
132  & xit(1)*gauss2d6(2,i)+
133  & xit(2)*(1-gauss2d6(1,i)-gauss2d6(2,i))
134 
135  etl= etat(3)*gauss2d6(1,i)+
136  & etat(1)*gauss2d6(2,i)+
137  & etat(2)*(1-gauss2d6(1,i)-gauss2d6(2,i))
138 !
139  call evalshapefunc(xil,etl,xns,nopes,xnl)
140  call evalshapefunc(xil,etl,xl2s,nopes,ps)
141 !
142  nintpoint=nintpoint+1
143 !
144 ! projection of the integration point in the mean
145 ! slave plane onto the slave surface
146 !
147 ! projection of the master integration point onto the
148 ! master surface in order to get the local coordinates
149 ! own xn for every integration point?
150 !
151  call attachline(xl2m,ps,nnodelem,xilm,etlm,xn)
152  call evalshapefunc(xilm,etlm,xl2m,nnodelem,pm)
153 !
154  pslavsurf(1,nintpoint)=xil
155  pslavsurf(2,nintpoint)=etl
156 !
157 ! weights sum up to 0.5 for triangles in local coordinates
158 !
159  pslavsurf(3,nintpoint)=2.d0*area*weight2d6(i)
160  enddo
161  enddo
162 !
163  return
subroutine sutherland_hodgman(nopes, xn, xl2sp, xl2mp, nodem, ipe, ime, iactiveline, nactiveline, ifreeintersec, nelemm, nnodelem, nvertex, pvertex)
Definition: sutherland_hodgman.f:26
subroutine attachline(pneigh, pnode, nterms, xil, etl, xn)
Definition: attachline.f:20
subroutine evalshapefunc(xil, etl, xl2, nopes, p)
Definition: evalshapefunc.f:20
integer function modf(m, ix)
Definition: modf.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)