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

Go to the source code of this file.

Functions/Subroutines

subroutine tiefaccont (lakon, ipkon, kon, ntie, tieset, nset, set, istartset, iendset, ialset, itiefac, islavsurf, islavnode, imastnode, nslavnode, nmastnode, nslavs, nmasts, ifacecount, iponoels, inoels, ifreenoels, mortar, ipoface, nodface, nk, xnoels)
 

Function/Subroutine Documentation

◆ tiefaccont()

subroutine tiefaccont ( character*8, dimension(*)  lakon,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
integer  ntie,
character*81, dimension(3,*)  tieset,
integer  nset,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer, dimension(2,*)  itiefac,
integer, dimension(2,*)  islavsurf,
integer, dimension(*)  islavnode,
integer, dimension(*)  imastnode,
integer, dimension(ntie+1)  nslavnode,
integer, dimension(ntie+1)  nmastnode,
integer  nslavs,
integer  nmasts,
integer  ifacecount,
integer, dimension(*)  iponoels,
integer, dimension(2,*)  inoels,
integer  ifreenoels,
integer  mortar,
integer, dimension(*)  ipoface,
integer, dimension(5,*)  nodface,
integer  nk,
real*8, dimension(*)  xnoels 
)
23 !
24 ! Catalogueing the slave faces (itieface, islavsurf)
25 ! the slave nodes (islavnode, nslavnode)
26 ! the slave faces to which the slave nodes
27 ! belong
28 ! the master nodes (imastnode, nmastnode; only
29 ! for surface-to-surface contact)
30 !
31 ! Authors: Li,Yang; Rakotonanahary, Samoela;
32 !
33  implicit none
34 !
35  logical nodeslavsurf
36 !
37  character*8 lakon(*)
38  character*81 tieset(3,*),slavset,mastset,set(*)
39 !
40  logical exist
41 !
42  integer ntie,i,j,k,l,nset,istartset(*),iendset(*),ialset(*),
43  & ifaces,nelems,jfaces,ifacem,nelemm,nslavs,nmasts,jface,
44  & jfacem,indexe,nopes,nopem,ipkon(*),kon(*),id,nodef(9),
45  & ifaceq(8,6),ifacet(6,4),ifacew1(4,5),ifacew2(8,5),node,
46  & itiefac(2,*),islavsurf(2,*),islavnode(*),imastnode(*),
47  & nslavnode(ntie+1),nmastnode(ntie+1),ifacecount,islav,imast,
48  & ipos,index1,iponoels(*),inoels(2,*),ifreenoels,ifreenoelold,
49  & mortar,numbern,numberf,iface,kflag,nk,ipoface(*),
50  & nodface(5,*),nface,nelem,nope
51 !
52  real*8 xnoels(*)
53 !
54 ! nslavnode: num of slave nodes
55 ! islavnode: all slave nodes, tie by tie, ordered within one tie constraint
56 ! nmastnode: num of master nodes
57 ! imastnode: all master nodes, tie by tie, ordered within one tie constraint
58 ! islavsurf: all slave faces
59 ! itiefac: pointer into field islavsurf
60 !
61 ! nodes per face for hex elements
62 !
63  data ifaceq /4,3,2,1,11,10,9,12,
64  & 5,6,7,8,13,14,15,16,
65  & 1,2,6,5,9,18,13,17,
66  & 2,3,7,6,10,19,14,18,
67  & 3,4,8,7,11,20,15,19,
68  & 4,1,5,8,12,17,16,20/
69 !
70 ! nodes per face for tet elements
71 !
72  data ifacet /1,3,2,7,6,5,
73  & 1,2,4,5,9,8,
74  & 2,3,4,6,10,9,
75  & 1,4,3,8,10,7/
76 !
77 ! nodes per face for linear wedge elements
78 !
79  data ifacew1 /1,3,2,0,
80  & 4,5,6,0,
81  & 1,2,5,4,
82  & 2,3,6,5,
83  & 3,1,4,6/
84 !
85 ! nodes per face for quadratic wedge elements
86 !
87  data ifacew2 /1,3,2,9,8,7,0,0,
88  & 4,5,6,10,11,12,0,0,
89  & 1,2,5,4,7,14,10,13,
90  & 2,3,6,5,8,15,11,14,
91  & 3,1,4,6,9,13,12,15/
92 !
93  ifacecount=0
94  nslavs=0
95  nmasts=0
96  ifreenoels=0
97  nodeslavsurf=.false.
98 !
99 ! counters for new fields islavsurf and itiefac
100 !
101  do i=1,ntie
102 !
103 ! check for contact conditions
104 !
105  if((tieset(1,i)(81:81).eq.'C').or.
106  & (tieset(1,i)(81:81).eq.'-')) then
107  slavset=tieset(2,i)
108 !
109 ! check whether facial slave surface;
110 !
111  ipos=index(slavset,' ')-1
112 !
113 ! default for node-to-surface contact is
114 ! a nodal slave surface
115 !
116  if(slavset(ipos:ipos).eq.'S') then
117  nodeslavsurf=.true.
118  endif
119 !
120 ! determining the slave surface
121 !
122  do j=1,nset
123  if(set(j).eq.slavset) exit
124  enddo
125  if(j.gt.nset) then
126  do j=1,nset
127  if((set(j)(1:ipos-1).eq.slavset(1:ipos-1)).and.
128  & (set(j)(ipos:ipos).eq.'T')) then
129  nodeslavsurf=.false.
130  exit
131  endif
132  enddo
133  endif
134 !
135  islav=j
136 !
137  if((mortar.eq.0).and.(nodeslavsurf)) then
138 !
139 ! nodal slave surface and node-to-surface contact
140 !
141 ! storing the slave nodes in islavnode (sorted)
142 !
143  nslavnode(i)=nslavs
144  numbern=0
145  do j=istartset(islav),iendset(islav)
146  if(ialset(j).gt.0) then
147  k=ialset(j)
148  call nident(islavnode(nslavs+1),k,numbern,id)
149  if(id.gt.0) then
150  if(islavnode(nslavs+id).eq.k) cycle
151  endif
152  numbern=numbern+1
153  do l=numbern,id+2,-1
154  islavnode(nslavs+l)=islavnode(nslavs+l-1)
155  enddo
156  islavnode(nslavs+id+1)=k
157  else
158  k=ialset(j-2)
159  do
160  k=k-ialset(j)
161  if(k.ge.ialset(j-1)) exit
162  call nident(islavnode(nslavs+1),k,numbern,id)
163  if(id.gt.0) then
164  if(islavnode(nslavs+id).eq.k) cycle
165  endif
166  numbern=numbern+1
167  do l=numbern,id+2,-1
168  islavnode(nslavs+l)=islavnode(nslavs+l-1)
169  enddo
170  islavnode(nslavs+id+1)=k
171  enddo
172  endif
173  nslavnode(i+1)=nslavnode(i)+numbern
174  enddo
175 !
176 ! check all external solid faces whether they contain
177 ! slave nodes
178 !
179 ! islavsurf(1,*) contains the faces (ordered)
180 ! islavsurf(2,*) contains the position of the
181 ! original order
182 !
183  itiefac(1,i)=ifacecount+1
184  numberf=0
185 !
186  do j=1,nk
187  index1=ipoface(j)
188  do
189  if(index1.eq.0) exit
190  iface=nodface(4,index1)
191  do k=0,3
192  if(k.eq.0) then
193  node=j
194  else
195  node=nodface(k,index1)
196  endif
197 !
198 ! check whether node belongs to slave surface
199 !
200  call nident(islavnode(nslavs+1),node,numbern,
201  & id)
202  if(id.gt.0) then
203  if(islavnode(nslavs+id).eq.node) then
204  call nident2(islavsurf(1,ifacecount+1),
205  & iface,numberf,id)
206 !
207  ipos=0
208  if(id.gt.0) then
209  if(islavsurf(1,ifacecount+id).eq.iface)
210  & then
211 ! original position (not ordered)
212  ipos=islavsurf(2,ifacecount+id)
213  endif
214  endif
215 !
216 ! check whether new face
217 !
218  if(ipos.eq.0) then
219  numberf=numberf+1
220  do l=ifacecount+numberf,ifacecount+id+2
221  & ,-1
222  islavsurf(1,l)=islavsurf(1,l-1)
223  islavsurf(2,l)=islavsurf(2,l-1)
224  enddo
225  islavsurf(1,ifacecount+id+1)=iface
226  islavsurf(2,ifacecount+id+1)=numberf
227 ! original position (not ordered)
228  ipos=numberf
229  endif
230 !
231 ! update info to which faces a slave node
232 ! belongs
233 !
234  ifreenoelold=iponoels(node)
235  ifreenoels=ifreenoels+1
236  iponoels(node)=ifreenoels
237  inoels(1,ifreenoels)=ifacecount+ipos
238 !
239 ! determining the area coefficient (= the
240 ! force coefficient corresponding to a
241 ! uniform pressure on the face)
242 !
243  nelem=int(iface/10.d0)
244  jface=iface-10*nelem
245 !
246  indexe=ipkon(nelem)
247  if(lakon(nelem)(4:4).eq.'2') then
248  nopes=8
249  nface=6
250  elseif(lakon(nelem)(4:4).eq.'8') then
251  nopes=4
252  nface=6
253  elseif(lakon(nelem)(4:5).eq.'10') then
254  nopes=6
255  nface=4
256  elseif(lakon(nelem)(4:4).eq.'4') then
257  nopes=3
258  nface=4
259  elseif(lakon(nelem)(4:5).eq.'15') then
260  if(jface.le.2) then
261  nopes=6
262  else
263  nopes=8
264  endif
265  nface=5
266  nope=15
267  elseif(lakon(nelem)(4:4).eq.'6') then
268  if(jface.le.2) then
269  nopes=3
270  else
271  nopes=4
272  endif
273  nface=5
274  nope=6
275  else
276  cycle
277  endif
278 !
279 ! determining the nodes of the face
280 !
281  if(nface.eq.4) then
282  do l=1,nopes
283  nodef(l)=kon(indexe+ifacet(l,jface))
284  if(nodef(l).eq.node) then
285  if(nopes.eq.3) then
286  xnoels(ifreenoels)=1.d0/3.d0
287 !
288 ! for 6-node faces the weights 0 and 1/3 are changed into
289 ! 1/999 and 332/999 in order to avoid problems with the
290 ! zero-area-check in springforc_n2f*f and springstiff_n2f*f
291 !
292  elseif(l.le.3) then
293  xnoels(ifreenoels)=1.d0/999.d0
294  else
295  xnoels(ifreenoels)
296  & =332.d0/999.d0
297  endif
298  endif
299  enddo
300  elseif(nface.eq.5) then
301  if(nope.eq.6) then
302  do l=1,nopes
303  nodef(l)=
304  & kon(indexe+ifacew1(l,jface))
305  if(nodef(l).eq.node) then
306  if(nopes.eq.3) then
307  xnoels(ifreenoels)=1.d0/3.d0
308  else
309  xnoels(ifreenoels)=1.d0/4.d0
310  endif
311  endif
312  enddo
313  elseif(nope.eq.15) then
314  do l=1,nopes
315  nodef(l)=
316  & kon(indexe+ifacew2(l,jface))
317  if(nodef(l).eq.node) then
318  if(nopes.eq.6) then
319  if(l.le.3) then
320  xnoels(ifreenoels)=
321  & 1.d0/999.d0
322  else
323  xnoels(ifreenoels)=
324  & 332.d0/999.d0
325  endif
326  else
327 !
328 ! for a 8-node face a distribution of 1/100 at the vertex
329 ! nodes and 24/100 at the midnodes (instead of -1/12 and 1/3)
330 !
331  if(l.le.4) then
332  xnoels(ifreenoels)=
333  & 1.d0/100.d0
334  else
335  xnoels(ifreenoels)=
336  & 24.d0/100.d0
337  endif
338  endif
339  endif
340  enddo
341  endif
342  elseif(nface.eq.6) then
343  do l=1,nopes
344  nodef(l)=kon(indexe+ifaceq(l,jface))
345  if(nodef(l).eq.node) then
346  if(nopes.eq.4) then
347  xnoels(ifreenoels)=1.d0/4.d0
348 !
349 ! for a 8-node face a distribution of 1/100 at the vertex
350 ! nodes and 24/100 at the midnodes (instead of -1/12 and 1/3)
351 !
352  elseif(l.le.4) then
353  xnoels(ifreenoels)=1.d0/100.d0
354  else
355  xnoels(ifreenoels)=
356  & 24.d0/100.d0
357  endif
358 
359  endif
360  enddo
361  endif
362 !
363  inoels(2,ifreenoels)=ifreenoelold
364  endif
365  endif
366  enddo
367  index1=nodface(5,index1)
368  enddo
369  enddo
370 !
371 ! restoring the right order of islavsurf(1,*);
372 ! thereafter islavsurf(2,*) is obsolete for node-to-surface
373 ! contact
374 !
375  kflag=2
376  call isorti(islavsurf(1,ifacecount+1),numberf,kflag)
377 !
378 ! update ifacecount and itiefac
379 !
380  itiefac(2,i)=ifacecount+numberf
381 !
382  nslavs=nslavnode(i+1)
383  ifacecount=itiefac(2,i)
384 !
385  else
386 !
387 ! element face slave surface (node-to-surface or
388 ! surface-to-surface contact)
389 !
390  nslavnode(i)=nslavs
391 !
392  itiefac(1,i)=ifacecount+1
393  do j=istartset(islav),iendset(islav)
394  if(ialset(j).gt.0) then
395 !
396 ! store the slave face in islavsurf
397 !
398  ifacecount=ifacecount+1
399  islavsurf(1,ifacecount)=ialset(j)
400 !
401 ! store the nodes belonging to the slave face
402 ! in islavnode
403 !
404  ifaces=ialset(j)
405  nelems=int(ifaces/10)
406  jfaces=ifaces - nelems*10
407  indexe=ipkon(nelems)
408 !
409  if(lakon(nelems)(4:5).eq.'20') then
410  nopes=8
411  elseif(lakon(nelems)(4:4).eq.'8') then
412  nopes=4
413  elseif(lakon(nelems)(4:5).eq.'10') then
414  nopes=6
415  elseif(lakon(nelems)(4:4).eq.'4') then
416  nopes=3
417  endif
418 !
419  if(lakon(nelems)(4:4).eq.'6') then
420  if(jfaces.le.2) then
421  nopes=3
422  else
423  nopes=4
424  endif
425  endif
426  if(lakon(nelems)(4:5).eq.'15') then
427  if(jfaces.le.2) then
428  nopes=6
429  else
430  nopes=8
431  endif
432  endif
433 !
434  do l=1,nopes
435  if((lakon(nelems)(4:4).eq.'2').or.
436  & (lakon(nelems)(4:4).eq.'8')) then
437  node=kon(indexe+ifaceq(l,jfaces))
438  elseif((lakon(nelems)(4:4).eq.'4').or.
439  & (lakon(nelems)(4:5).eq.'10')) then
440  node=kon(indexe+ifacet(l,jfaces))
441  elseif(lakon(nelems)(4:4).eq.'6') then
442  node=kon(indexe+ifacew1(l,jfaces))
443  elseif(lakon(nelems)(4:5).eq.'15') then
444  node=kon(indexe+ifacew2(l,jfaces))
445  endif
446  call nident(islavnode(nslavnode(i)+1),node,
447  & nslavs-nslavnode(i),id)
448  exist=.false.
449  if(id.gt.0) then
450  if(islavnode(nslavnode(i)+id).eq.node) then
451  exist=.true.
452  endif
453  endif
454  if(.not.exist) then
455  nslavs=nslavs+1
456  do k=nslavs,nslavnode(i)+id+2,-1
457  islavnode(k)=islavnode(k-1)
458  enddo
459  islavnode(nslavnode(i)+id+1)=node
460  endif
461 !
462 ! filling fields iponoels and inoels
463 !
464  ifreenoelold=iponoels(node)
465  ifreenoels=ifreenoels+1
466  iponoels(node)=ifreenoels
467  inoels(1,ifreenoels)=ifacecount
468 !
469 ! filling xnoels with the coefficient corresponding
470 ! to a constant pressure (sum over all nodes must
471 ! be 1)
472 !
473  if(nopes.eq.3) then
474  xnoels(ifreenoels)=1.d0/3.d0
475  elseif(nopes.eq.4) then
476  xnoels(ifreenoels)=1.d0/4.d0
477  elseif(nopes.eq.6) then
478  if(l.le.3) then
479  xnoels(ifreenoels)=1.d0/999.d0
480  else
481  xnoels(ifreenoels)=332.d0/999.d0
482  endif
483  elseif(nopes.eq.7) then
484  if(l.le.3) then
485  xnoels(ifreenoels)=1.d0/20.d0
486  elseif(l.le.6) then
487  xnoels(ifreenoels)=2.d0/15.d0
488  else
489  xnoels(ifreenoels)=9.d0/20.d0
490  endif
491  elseif(nopes.eq.8) then
492 !
493 ! for a 8-node face a distribution of 1/100 at the vertex
494 ! nodes and 24/100 at the midnodes (instead of -1/12 and 1/3)
495 !
496  if(l.le.4) then
497  xnoels(ifreenoels)=1.d0/100.d0
498  else
499  xnoels(ifreenoels)=24.d0/100.d0
500  endif
501  elseif(nopes.eq.9) then
502  if(l.le.4) then
503  xnoels(ifreenoels)=1.d0/36.d0
504  elseif(l.le.8) then
505  xnoels(ifreenoels)=1.d0/9.d0
506  else
507  xnoels(ifreenoels)=4.d0/9.d0
508  endif
509  endif
510  inoels(2,ifreenoels)=ifreenoelold
511  enddo
512 !
513  endif
514  enddo
515  nslavnode(i+1)=nslavs
516  itiefac(2,i)=ifacecount
517  endif
518 !
519 ! determining the master surface
520 !
521  mastset=tieset(3,i)
522  do j=1,nset
523  if(set(j).eq.mastset) exit
524  enddo
525  if(j.gt.nset) then
526  write(*,*) '*ERROR in tiefaccont: master surface'
527  write(*,*) ' does not exist'
528  call exit(201)
529  endif
530  imast=j
531  nmastnode(i)=nmasts
532 !
533  do j=istartset(imast),iendset(imast)
534 !
535 ! create imastnode, and nmastnode
536 !
537  ifacem=ialset(j)
538  nelemm=int(ifacem/10)
539  jfacem=ifacem - nelemm*10
540  indexe=ipkon(nelemm)
541 !
542  if(lakon(nelemm)(4:5).eq.'20') then
543  nopem=8
544  elseif(lakon(nelemm)(4:4).eq.'8') then
545  nopem=4
546  elseif(lakon(nelemm)(4:5).eq.'10') then
547  nopem=6
548  elseif(lakon(nelemm)(4:4).eq.'4') then
549  nopem=3
550  endif
551 !
552  if(lakon(nelemm)(4:4).eq.'6') then
553  if(jfacem.le.2) then
554  nopem=3
555  else
556  nopem=4
557  endif
558  endif
559  if(lakon(nelemm)(4:5).eq.'15') then
560  if(jfacem.le.2) then
561  nopem=6
562  else
563  nopem=8
564  endif
565  endif
566 !
567  do l=1,nopem
568  if((lakon(nelemm)(4:4).eq.'2').or.
569  & (lakon(nelemm)(4:4).eq.'8')) then
570  node=kon(indexe+ifaceq(l,jfacem))
571  elseif((lakon(nelemm)(4:4).eq.'4').or.
572  & (lakon(nelemm)(4:5).eq.'10')) then
573  node=kon(indexe+ifacet(l,jfacem))
574  elseif(lakon(nelemm)(4:4).eq.'6') then
575  node=kon(indexe+ifacew1(l,jfacem))
576  elseif(lakon(nelemm)(4:5).eq.'15') then
577  node=kon(indexe+ifacew2(l,jfacem))
578  endif
579  call nident(imastnode(nmastnode(i)+1),node,
580  & nmasts-nmastnode(i),id)
581  exist=.false.
582  if(id.gt.0) then
583  if(imastnode(nmastnode(i)+id).eq.node) then
584  exist=.true.
585  endif
586  endif
587  if(exist) cycle
588  nmasts=nmasts+1
589  do k=nmasts,nmastnode(i)+id+2,-1
590  imastnode(k)=imastnode(k-1)
591  enddo
592  imastnode(nmastnode(i)+id+1)=node
593  enddo
594 !
595  enddo
596  nmastnode(i+1)=nmasts
597 !
598  else
599 !
600 ! no contact tie
601 !
602  nslavnode(i+1)=nslavnode(i)
603 c if(mortar.eq.1) nmastnode(i+1)=nmastnode(i)
604  nmastnode(i+1)=nmastnode(i)
605  endif
606  enddo
607 !
608  return
subroutine nident2(x, px, n, id)
Definition: nident2.f:27
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine isorti(ix, n, kflag)
Definition: isorti.f:6
Hosted by OpenAircraft.com, (Michigan UAV, LLC)