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

Go to the source code of this file.

Functions/Subroutines

subroutine adjustcontactnodes (tieset, ntie, itietri, cg, straight, co, vold, xo, yo, zo, x, y, z, nx, ny, nz, istep, iinc, iit, mi, imastop, nslavnode, islavnode, set, nset, istartset, iendset, ialset, tietol, clearini, clearslavnode, itiefac, ipkon, kon, lakon, islavsurf)
 

Function/Subroutine Documentation

◆ adjustcontactnodes()

subroutine adjustcontactnodes ( character*81, dimension(3,*)  tieset,
integer  ntie,
integer, dimension(2,ntie)  itietri,
real*8, dimension(3,*)  cg,
real*8, dimension(16,*)  straight,
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  istep,
integer  iinc,
integer  iit,
integer, dimension(*)  mi,
integer, dimension(3,*)  imastop,
integer, dimension(*)  nslavnode,
integer, dimension(*)  islavnode,
character*81, dimension(*)  set,
integer  nset,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
real*8, dimension(3,*)  tietol,
real*8, dimension(3,9,*)  clearini,
real*8, dimension(3,*)  clearslavnode,
integer, dimension(2,*)  itiefac,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
integer, dimension(2,*)  islavsurf 
)
24 !
25 ! Adjusting contact nodes if requested by the user
26 !
27  implicit none
28 !
29  character*8 lakon(*)
30  character*81 tieset(3,*),slavset,set(*),noset
31 !
32  integer ntie,
33  & itietri(2,ntie),node,neigh(1),kneigh,i,j,k,l,m,isol,iset,
34  & idummy,itri,ll,kflag,n,nx(*),ny(*),istep,iinc,mi(*),
35  & nz(*),nstart,iit,imastop(3,*),itriangle(100),ntriangle,
36  & ntriangle_,itriold,itrinew,id,nslavnode(*),islavnode(*),
37  & ipos,nset,istartset(*),iendset(*),ialset(*),iclear,
38  & istart,ilength,nope,nopes,nelems,jj,ifaces,itiefac(2,*),
39  & jfaces,islavsurf(2,*),ifaceq(8,6),ifacet(6,4),ifacew1(4,5),
40  & ifacew2(8,5),ipkon(*),kon(*)
41 !
42  real*8 cg(3,*),straight(16,*),co(3,*),vold(0:mi(2),*),p(3),
43  & dist,xo(*),yo(*),zo(*),x(*),y(*),z(*),tietol(3,*),adjust,
44  & clearini(3,9,*),clearslavnode(3,*),clear
45 !
46 ! nodes per face for hex elements
47 !
48  data ifaceq /4,3,2,1,11,10,9,12,
49  & 5,6,7,8,13,14,15,16,
50  & 1,2,6,5,9,18,13,17,
51  & 2,3,7,6,10,19,14,18,
52  & 3,4,8,7,11,20,15,19,
53  & 4,1,5,8,12,17,16,20/
54 !
55 ! nodes per face for tet elements
56 !
57  data ifacet /1,3,2,7,6,5,
58  & 1,2,4,5,9,8,
59  & 2,3,4,6,10,9,
60  & 1,4,3,8,10,7/
61 !
62 ! nodes per face for linear wedge elements
63 !
64  data ifacew1 /1,3,2,0,
65  & 4,5,6,0,
66  & 1,2,5,4,
67  & 2,3,6,5,
68  & 3,1,4,6/
69 !
70 ! nodes per face for quadratic wedge elements
71 !
72  data ifacew2 /1,3,2,9,8,7,0,0,
73  & 4,5,6,10,11,12,0,0,
74  & 1,2,5,4,7,14,10,13,
75  & 2,3,6,5,8,15,11,14,
76  & 3,1,4,6,9,13,12,15/
77 !
78  do i=1,ntie
79  if(tieset(1,i)(81:81).ne.'C') cycle
80  kneigh=1
81  slavset=tieset(2,i)
82 !
83 ! check whether a clearance has been specified by the user
84 !
85 c tietol(3,i)=1.2357111317d0
86  if((tietol(3,i).gt.1.2357111316d0).and.
87  & (tietol(3,i).lt.1.2357111318d0)) then
88  iclear=0
89  else
90  iclear=1
91  clear=tietol(3,i)
92  endif
93 !
94 ! check whether an adjust node set has been defined
95 ! only checked in the first increment of the first step
96 !
97  if(istep.eq.1) then
98  iset=0
99  if(tieset(1,i)(1:1).ne.' ') then
100  noset(1:80)=tieset(1,i)(1:80)
101  noset(81:81)=' '
102  ipos=index(noset,' ')
103  noset(ipos:ipos)='N'
104  do iset=1,nset
105  if(set(iset).eq.noset) exit
106  enddo
107  kflag=1
108  call isortii(ialset(istartset(iset)),idummy,
109  & iendset(iset)-istartset(iset)+1,kflag)
110  endif
111  endif
112 !
113 ! search a master face for each slave node; determine the
114 ! distance and adjust it if requested by the user
115 !
116  nstart=itietri(1,i)-1
117  n=itietri(2,i)-nstart
118  if(n.lt.kneigh) kneigh=n
119  do j=1,n
120  xo(j)=cg(1,nstart+j)
121  x(j)=xo(j)
122  nx(j)=j
123  yo(j)=cg(2,nstart+j)
124  y(j)=yo(j)
125  ny(j)=j
126  zo(j)=cg(3,nstart+j)
127  z(j)=zo(j)
128  nz(j)=j
129  enddo
130  kflag=2
131  call dsort(x,nx,n,kflag)
132  call dsort(y,ny,n,kflag)
133  call dsort(z,nz,n,kflag)
134 !
135  do j=nslavnode(i)+1,nslavnode(i+1)
136  node=islavnode(j)
137 !
138  do k=1,3
139  p(k)=co(k,node)
140 c p(k)=co(k,node)+vold(k,node)
141  enddo
142 !
143 ! determining the kneigh neighboring master contact
144 ! triangle centers of gravity
145 !
146  call near3d(xo,yo,zo,x,y,z,nx,ny,nz,p(1),p(2),p(3),
147  & n,neigh,kneigh)
148 !
149  isol=0
150 !
151  itriold=0
152  itri=neigh(1)+itietri(1,i)-1
153  ntriangle=0
154  ntriangle_=100
155 !
156  loop1: do
157  do l=1,3
158  ll=4*l-3
159  dist=straight(ll,itri)*p(1)+
160  & straight(ll+1,itri)*p(2)+
161  & straight(ll+2,itri)*p(3)+
162  & straight(ll+3,itri)
163  if(dist.gt.1.d-6) then
164  itrinew=imastop(l,itri)
165  if(itrinew.eq.0) then
166 c write(*,*) '**border reached'
167  exit loop1
168  elseif(itrinew.eq.itriold) then
169 c write(*,*) '**solution in between triangles'
170  isol=itri
171  exit loop1
172  else
173  call nident(itriangle,itrinew,ntriangle,id)
174  if(id.gt.0) then
175  if(itriangle(id).eq.itrinew) then
176 c write(*,*) '**circular path; no solution'
177  exit loop1
178  endif
179  endif
180  ntriangle=ntriangle+1
181  if(ntriangle.gt.ntriangle_) then
182 c write(*,*) '**too many iterations'
183  exit loop1
184  endif
185  do k=ntriangle,id+2,-1
186  itriangle(k)=itriangle(k-1)
187  enddo
188  itriangle(id+1)=itrinew
189  itriold=itri
190  itri=itrinew
191  cycle loop1
192  endif
193  elseif(l.eq.3) then
194 c write(*,*) '**regular solution'
195  isol=itri
196  exit loop1
197  endif
198  enddo
199  enddo loop1
200 !
201 ! calculating the distance
202 !
203  if(isol.ne.0) then
204  dist=straight(13,itri)*p(1)+
205  & straight(14,itri)*p(2)+
206  & straight(15,itri)*p(3)+
207  & straight(16,itri)
208 !
209 ! check for an adjust parameter (only in the first
210 ! increment of the first step)
211 !
212  if(istep.eq.1) then
213  if(iset.ne.0) then
214 !
215 ! check whether node belongs to the adjust node
216 ! set
217 !
218  call nident(ialset(istartset(iset)),node,
219  & iendset(iset)-istartset(iset)+1,id)
220  if(id.gt.0) then
221  if(ialset(istartset(iset)+id-1).eq.node) then
222  do k=1,3
223  co(k,node)=co(k,node)-
224  & dist*straight(12+k,itri)
225  enddo
226  dist=0.d0
227  endif
228  endif
229  elseif(dabs(tietol(1,i)).ge.2.d0) then
230 !
231 ! adjust parameter
232 !
233  adjust=dabs(tietol(1,i))-2.d0
234 
235  if(dist.le.adjust) then
236  do k=1,3
237  co(k,node)=co(k,node)-
238  & dist*straight(12+k,itri)
239  enddo
240  dist=0.d0
241  endif
242  endif
243  endif
244 !
245 ! clearance in each slave node
246 !
247  if(iclear.eq.1) then
248  do k=1,3
249  clearslavnode(k,j)=(clear-dist)*straight(12+k,itri)
250  enddo
251  endif
252 !
253  endif
254  enddo
255 !
256 ! loop over all slave faces
257 !
258  if(iclear.eq.1) then
259  do jj=itiefac(1,i), itiefac(2,i)
260  ifaces=islavsurf(1,jj)
261  nelems=int(ifaces/10)
262  jfaces=ifaces - nelems*10
263 !
264  if(lakon(nelems)(4:5).eq.'8R') then
265  nopes=4
266  nope=8
267  elseif(lakon(nelems)(4:4).eq.'8') then
268  nopes=4
269  nope=8
270  elseif(lakon(nelems)(4:6).eq.'20R') then
271  nopes=8
272  nope=20
273  elseif(lakon(nelems)(4:5).eq.'20') then
274  nopes=8
275  nope=20
276  elseif(lakon(nelems)(4:5).eq.'10') then
277  nopes=6
278  nope=10
279  elseif(lakon(nelems)(4:4).eq.'4') then
280  nopes=3
281  nope=4
282 !
283 ! treatment of wedge faces
284 !
285  elseif(lakon(nelems)(4:4).eq.'6') then
286  nope=6
287  if(jfaces.le.2) then
288  nopes=3
289  else
290  nopes=4
291  endif
292  elseif(lakon(nelems)(4:5).eq.'15') then
293  nope=15
294  if(jfaces.le.2) then
295  nopes=6
296  else
297  nopes=8
298  endif
299  endif
300 !
301 ! actual position of the nodes belonging to the
302 ! slave surface
303 !
304  istart=nslavnode(i)+1
305  ilength=nslavnode(i+1)-nslavnode(i)
306 !
307  do m=1,nopes
308  if((nope.eq.20).or.(nope.eq.8))then
309  node=kon(ipkon(nelems)+ifaceq(m,jfaces))
310  elseif((nope.eq.10).or.(nope.eq.4).or.(nope.eq.14))
311  & then
312  node=kon(ipkon(nelems)+ifacet(m,jfaces))
313  elseif(nope.eq.15) then
314  node=kon(ipkon(nelems)+ifacew2(m,jfaces))
315  else
316  node=kon(ipkon(nelems)+ifacew1(m,jfaces))
317  endif
318  call nident(islavnode(istart),node,ilength,id)
319  do k=1,3
320  clearini(k,m,jj)=
321  & clearslavnode(k,nslavnode(i)+id)
322  enddo
323  enddo
324  enddo
325  endif
326 !
327 ! transfer the clearance from the slave nodes in islavnode
328 ! to the nodes per slave face in islavsurf
329 !
330  enddo
331 !
332  return
subroutine near3d(xo, yo, zo, x, y, z, nx, ny, nz, xp, yp, zp, n, neighbor, k)
Definition: near3d.f:20
subroutine isortii(ix, iy, n, kflag)
Definition: isortii.f:6
static double * dist
Definition: radflowload.c:42
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine dsort(dx, iy, n, kflag)
Definition: dsort.f:6
Hosted by OpenAircraft.com, (Michigan UAV, LLC)