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

Go to the source code of this file.

Functions/Subroutines

subroutine gencontelem_n2f (tieset, ntie, itietri, ne, ipkon, kon, lakon, cg, straight, ifree, koncont, co, vold, xo, yo, zo, x, y, z, nx, ny, nz, ielmat, elcon, istep, iinc, iit, ncmat_, ntmat_, nmethod, mi, imastop, nslavnode, islavnode, islavsurf, itiefac, areaslav, iponoels, inoels, springarea, set, nset, istartset, iendset, ialset, tietol, reltime, filab, nasym, xnoels, icutb, ne0, jobnamef)
 

Function/Subroutine Documentation

◆ gencontelem_n2f()

subroutine gencontelem_n2f ( character*81, dimension(3,*)  tieset,
integer  ntie,
integer, dimension(2,ntie)  itietri,
integer  ne,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
real*8, dimension(3,*)  cg,
real*8, dimension(16,*)  straight,
integer  ifree,
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(mi(3),*)  ielmat,
real*8, dimension(0:ncmat_,ntmat_,*)  elcon,
integer  istep,
integer  iinc,
integer  iit,
integer  ncmat_,
integer  ntmat_,
integer  nmethod,
integer, dimension(*)  mi,
integer, dimension(3,*)  imastop,
integer, dimension(*)  nslavnode,
integer, dimension(*)  islavnode,
integer, dimension(2,*)  islavsurf,
integer, dimension(2,*)  itiefac,
real*8, dimension(*)  areaslav,
integer, dimension(*)  iponoels,
integer, dimension(2,*)  inoels,
real*8, dimension(2,*)  springarea,
character*81, dimension(*)  set,
integer  nset,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
real*8, dimension(3,*)  tietol,
real*8  reltime,
character*87, dimension(*)  filab,
integer  nasym,
real*8, dimension(*)  xnoels,
integer  icutb,
integer  ne0,
character*132, dimension(*)  jobnamef 
)
27 !
28 ! generate contact elements for the slave contact nodes
29 !
30  implicit none
31 !
32  character*8 lakon(*)
33  character*33 cfile
34  character*81 tieset(3,*),slavset,set(*),noset,setname
35  character*87 filab(*)
36  character*132 jobnamef(*)
37 !
38  integer ntie,ifree,nasym,icutb,ne0,
39  & itietri(2,ntie),ipkon(*),kon(*),koncont(4,*),ne,node,
40  & neigh(1),iflag,kneigh,i,j,k,l,isol,iset,idummy,
41  & itri,ll,kflag,n,nx(*),ny(*),istep,iinc,mi(*),
42  & nz(*),nstart,ielmat(mi(3),*),imat,ifaceq(8,6),ifacet(6,4),
43  & ifacew1(4,5),ifacew2(8,5),nelem,jface,indexe,iit,
44  & nface,nope,nodef(9),ncmat_,ntmat_,index1,indexel,
45  & nmethod,iteller,ifaces,jfaces,irelslavface,number(4),lenset,
46  & imastop(3,*), itriangle(100),ntriangle,ntriangle_,itriold,
47  & itrinew,id,nslavnode(*),islavnode(*),islavsurf(2,*),
48  & itiefac(2,*),iponoels(*),inoels(2,*),konl(26),nelems,m,
49  & mint2d,nopes,ipos,nset,istartset(*),iendset(*),
50  & ialset(*)
51 !
52  real*8 cg(3,*),straight(16,*),co(3,*),vold(0:mi(2),*),p(3),
53  & dist,xo(*),yo(*),zo(*),x(*),y(*),z(*),c0coef,
54  & beta,c0,elcon(0:ncmat_,ntmat_,*),weight,
55  & areaslav(*),springarea(2,*),xl2(3,9),area,xi,et,shp2(7,9),
56  & xs2(3,2),xsj2(3),adjust,tietol(3,*),reltime,
57  & clear,ratio(9),pl(3,9),
58  & pproj(3),al(3),xn(3),xm(3),dm,xnoels(*)
59 !
60 ! nodes per face for hex elements
61 !
62  data ifaceq /4,3,2,1,11,10,9,12,
63  & 5,6,7,8,13,14,15,16,
64  & 1,2,6,5,9,18,13,17,
65  & 2,3,7,6,10,19,14,18,
66  & 3,4,8,7,11,20,15,19,
67  & 4,1,5,8,12,17,16,20/
68 !
69 ! nodes per face for tet elements
70 !
71  data ifacet /1,3,2,7,6,5,
72  & 1,2,4,5,9,8,
73  & 2,3,4,6,10,9,
74  & 1,4,3,8,10,7/
75 !
76 ! nodes per face for linear wedge elements
77 !
78  data ifacew1 /1,3,2,0,
79  & 4,5,6,0,
80  & 1,2,5,4,
81  & 2,3,6,5,
82  & 3,1,4,6/
83 !
84 ! nodes per face for quadratic wedge elements
85 !
86  data ifacew2 /1,3,2,9,8,7,0,0,
87  & 4,5,6,10,11,12,0,0,
88  & 1,2,5,4,7,14,10,13,
89  & 2,3,6,5,8,15,11,14,
90  & 3,1,4,6,9,13,12,15/
91 !
92 ! flag for shape functions
93 !
94  data iflag /2/
95  data indexel /0/
96 !
97  save indexel
98 !
99  include "gauss.f"
100 !
101 ! opening a file to store the contact spring elements
102 !
103  if(filab(1)(3:3).eq.'C') then
104  iteller=iteller+1
105 !
106  do i=1,132
107  if(jobnamef(1)(i:i).eq.' ') exit
108  enddo
109  i=i-1
110  cfile=jobnamef(1)(1:i)//'.cel'
111  open(27,file=cfile,status='unknown',position='append')
112 !
113  setname(1:15)='contactelements'
114  lenset=15
115  number(1)=istep
116  number(2)=iinc
117  number(3)=icutb+1
118  number(4)=iit
119  if(number(4).le.0) number(4)=1
120  do i=1,4
121  setname(lenset+1:lenset+1)='_'
122  if(i.eq.1) then
123  setname(lenset+2:lenset+3)='st'
124  elseif(i.eq.2) then
125  setname(lenset+2:lenset+3)='in'
126  elseif(i.eq.3) then
127  setname(lenset+2:lenset+3)='at'
128  else
129  setname(lenset+2:lenset+3)='it'
130  endif
131  if(number(i).lt.10) then
132  write(setname(lenset+4:lenset+4),'(i1)') number(i)
133  lenset=lenset+4
134  elseif(number(i).lt.100) then
135  write(setname(lenset+4:lenset+5),'(i2)') number(i)
136  lenset=lenset+5
137  elseif(number(i).lt.1000) then
138  write(setname(lenset+4:lenset+6),'(i3)') number(i)
139  lenset=lenset+6
140  else
141  write(*,*) '*ERROR in gencontelem_f2f: no more than 1000'
142  write(*,*) ' steps/increments/cutbacks/iterations'
143  write(*,*) ' allowed (for output in '
144  write(*,*) ' contactelements.inp)'
145  stop
146  endif
147  enddo
148  endif
149 !
150  do i=1,ntie
151  if(tieset(1,i)(81:81).ne.'C') cycle
152  kneigh=1
153  slavset=tieset(2,i)
154  imat=int(tietol(2,i))
155  c0coef=elcon(4,1,imat)
156 !
157 ! check whether an adjust node set has been defined
158 ! only checked at the start of the first step
159 !
160 c if((istep.eq.1).and.(iinc.eq.1).and.(iit.le.0)) then
161  if((istep.eq.1).and.(iit.lt.0)) then
162  iset=0
163  if(tieset(1,i)(1:1).ne.' ') then
164  noset(1:80)=tieset(1,i)(1:80)
165  noset(81:81)=' '
166  ipos=index(noset,' ')
167  noset(ipos:ipos)='N'
168  do iset=1,nset
169  if(set(iset).eq.noset) exit
170  enddo
171  kflag=1
172  call isortii(ialset(istartset(iset)),idummy,
173  & iendset(iset)-istartset(iset)+1,kflag)
174  endif
175  endif
176 !
177 ! determine the area of the slave surfaces
178 !
179  do l = itiefac(1,i), itiefac(2,i)
180  ifaces = islavsurf(1,l)
181  nelems = int(ifaces/10)
182  jfaces = ifaces - nelems*10
183 !
184 ! Decide on the max integration points number, just consider 2D situation
185 !
186  if(lakon(nelems)(4:5).eq.'8R') then
187  mint2d=1
188  nopes=4
189  nope=8
190  elseif(lakon(nelems)(4:4).eq.'8') then
191  mint2d=4
192  nopes=4
193  nope=8
194  elseif(lakon(nelems)(4:6).eq.'20R') then
195  mint2d=4
196  nopes=8
197  nope=20
198  elseif(lakon(nelems)(4:5).eq.'20') then
199  mint2d=9
200  nopes=8
201  nope=20
202  elseif(lakon(nelems)(4:5).eq.'10') then
203  mint2d=3
204  nopes=6
205  nope=10
206  elseif(lakon(nelems)(4:4).eq.'4') then
207  mint2d=1
208  nopes=3
209  nope=4
210 !
211 ! treatment of wedge faces
212 !
213  elseif(lakon(nelems)(4:4).eq.'6') then
214  mint2d=1
215  nope=6
216  if(jfaces.le.2) then
217  nopes=3
218  else
219  nopes=4
220  endif
221  elseif(lakon(nelems)(4:5).eq.'15') then
222  nope=15
223  if(jfaces.le.2) then
224  mint2d=3
225  nopes=6
226  else
227  mint2d=4
228  nopes=8
229  endif
230  endif
231 !
232 ! actual position of the nodes belonging to the
233 ! slave surface
234 !
235  do j=1,nope
236  konl(j)=kon(ipkon(nelems)+j)
237  enddo
238 !
239  if((nope.eq.20).or.(nope.eq.8)) then
240  do m=1,nopes
241  do j=1,3
242  xl2(j,m)=co(j,konl(ifaceq(m,jfaces)))+
243  & vold(j,konl(ifaceq(m,jfaces)))
244  enddo
245  enddo
246  elseif((nope.eq.10).or.(nope.eq.4)) then
247  do m=1,nopes
248  do j=1,3
249  xl2(j,m)=co(j,konl(ifacet(m,jfaces)))+
250  & vold(j,konl(ifacet(m,jfaces)))
251  enddo
252  enddo
253  elseif(nope.eq.15) then
254  do m=1,nopes
255  do j=1,3
256  xl2(j,m)=co(j,konl(ifacew2(m,jfaces)))+
257  & vold(j,konl(ifacew2(m,jfaces)))
258  enddo
259  enddo
260  else
261  do m=1,nopes
262  do j=1,3
263  xl2(j,m)=co(j,konl(ifacew1(m,jfaces)))+
264  & vold(j,konl(ifacew1(m,jfaces)))
265  enddo
266  enddo
267  endif
268 !
269 ! calculating the area of the slave face
270 !
271  area=0.d0
272  do m=1,mint2d
273  if((lakon(nelems)(4:5).eq.'8R').or.
274  & ((lakon(nelems)(4:4).eq.'6').and.(nopes.eq.4))) then
275  xi=gauss2d1(1,m)
276  et=gauss2d1(2,m)
277  weight=weight2d1(m)
278  elseif((lakon(nelems)(4:4).eq.'8').or.
279  & (lakon(nelems)(4:6).eq.'20R').or.
280  & ((lakon(nelems)(4:5).eq.'15').and.
281  & (nopes.eq.8))) then
282  xi=gauss2d2(1,m)
283  et=gauss2d2(2,m)
284  weight=weight2d2(m)
285  elseif(lakon(nelems)(4:4).eq.'2') then
286  xi=gauss2d3(1,m)
287  et=gauss2d3(2,m)
288  weight=weight2d3(m)
289  elseif((lakon(nelems)(4:5).eq.'10').or.
290  & ((lakon(nelems)(4:5).eq.'15').and.
291  & (nopes.eq.6))) then
292  xi=gauss2d5(1,m)
293  et=gauss2d5(2,m)
294  weight=weight2d5(m)
295  elseif((lakon(nelems)(4:4).eq.'4').or.
296  & ((lakon(nelems)(4:4).eq.'6').and.
297  & (nopes.eq.3))) then
298  xi=gauss2d4(1,m)
299  et=gauss2d4(2,m)
300  weight=weight2d4(m)
301  endif
302 !
303  if(nopes.eq.9) then
304  call shape9q(xi,et,xl2,xsj2,xs2,shp2,iflag)
305  elseif(nopes.eq.8) then
306  call shape8q(xi,et,xl2,xsj2,xs2,shp2,iflag)
307  elseif(nopes.eq.4) then
308  call shape4q(xi,et,xl2,xsj2,xs2,shp2,iflag)
309  elseif(nopes.eq.6) then
310  call shape6tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
311  elseif(nopes.eq.7) then
312  call shape7tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
313  else
314  call shape3tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
315  endif
316  area=area+weight*dsqrt(xsj2(1)**2+xsj2(2)**2+
317  & xsj2(3)**2)
318  enddo
319  areaslav(l)=area
320  enddo
321 !
322 ! search a master face for each slave node and generate a contact
323 ! spring element if successful
324 !
325  nstart=itietri(1,i)-1
326  n=itietri(2,i)-nstart
327  if(n.lt.kneigh) kneigh=n
328  do j=1,n
329  xo(j)=cg(1,nstart+j)
330  x(j)=xo(j)
331  nx(j)=j
332  yo(j)=cg(2,nstart+j)
333  y(j)=yo(j)
334  ny(j)=j
335  zo(j)=cg(3,nstart+j)
336  z(j)=zo(j)
337  nz(j)=j
338  enddo
339  kflag=2
340  call dsort(x,nx,n,kflag)
341  call dsort(y,ny,n,kflag)
342  call dsort(z,nz,n,kflag)
343 !
344  do j=nslavnode(i)+1,nslavnode(i+1)
345  if(iit.le.0) springarea(2,j)=0.d0
346  node=islavnode(j)
347 !
348 ! calculating the area corresponding to the
349 ! slave node; is made up of the area
350 ! of the neighboring slave faces
351 !
352  area=0.d0
353  index1=iponoels(node)
354  do
355  if(index1.eq.0) exit
356  irelslavface=inoels(1,index1)
357  if((itiefac(1,i).le.irelslavface).and.
358  & (irelslavface.le.itiefac(2,i))) then
359  area=area+areaslav(irelslavface)*
360  & xnoels(index1)
361  endif
362  index1=inoels(2,index1)
363  enddo
364 !
365  do k=1,3
366  p(k)=co(k,node)+vold(k,node)
367  enddo
368 !
369 ! determining the kneigh neighboring master contact
370 ! triangle centers of gravity
371 !
372  call near3d(xo,yo,zo,x,y,z,nx,ny,nz,p(1),p(2),p(3),
373  & n,neigh,kneigh)
374 !
375  isol=0
376 !
377  itriold=0
378  itri=neigh(1)+itietri(1,i)-1
379  ntriangle=0
380  ntriangle_=100
381 !
382  loop1: do
383  do l=1,3
384  ll=4*l-3
385  dist=straight(ll,itri)*p(1)+
386  & straight(ll+1,itri)*p(2)+
387  & straight(ll+2,itri)*p(3)+
388  & straight(ll+3,itri)
389 !
390 ! 1.d-6 was increased to 1.d-3 on 19/04/2012
391 ! this is important for 2d-calculations or
392 ! calculations for which structures fit exactly
393 ! at their boundaries
394 !
395  if(dist.gt.1.d-3*dsqrt(area)) then
396  itrinew=imastop(l,itri)
397  if(itrinew.eq.0) then
398 c write(*,*) '**border reached'
399  exit loop1
400  elseif((itrinew.lt.itietri(1,i)).or.
401  & (itrinew.gt.itietri(2,i))) then
402 c write(*,*) '**border reached'
403  exit loop1
404  elseif(itrinew.eq.itriold) then
405 c write(*,*) '**solution in between triangles'
406  isol=itri
407  exit loop1
408  else
409  call nident(itriangle,itrinew,ntriangle,id)
410  if(id.gt.0) then
411  if(itriangle(id).eq.itrinew) then
412 c write(*,*) '**circular path;no solution'
413  exit loop1
414  endif
415  endif
416  ntriangle=ntriangle+1
417  if(ntriangle.gt.ntriangle_) then
418 c write(*,*) '**too many iterations'
419  exit loop1
420  endif
421  do k=ntriangle,id+2,-1
422  itriangle(k)=itriangle(k-1)
423  enddo
424  itriangle(id+1)=itrinew
425  itriold=itri
426  itri=itrinew
427  cycle loop1
428  endif
429  elseif(l.eq.3) then
430 c write(*,*) '**regular solution'
431  isol=itri
432  exit loop1
433  endif
434  enddo
435  enddo loop1
436 !
437 ! check whether distance is larger than c0:
438 ! no element is generated
439 !
440  if(isol.ne.0) then
441 !
442 ! determining the clearance
443 !
444 ! identifying the element face to which the
445 ! triangle belongs
446 !
447  nelem=int(koncont(4,itri)/10.d0)
448  jface=koncont(4,itri)-10*nelem
449 !
450  indexe=ipkon(nelem)
451  if(lakon(nelem)(4:5).eq.'20') then
452  nopes=8
453  nface=6
454  elseif(lakon(nelem)(4:4).eq.'8') then
455  nopes=4
456  nface=6
457  elseif(lakon(nelem)(4:5).eq.'10') then
458  nopes=6
459  nface=4
460  elseif(lakon(nelem)(4:4).eq.'4') then
461  nopes=3
462  nface=4
463  elseif(lakon(nelem)(4:5).eq.'15') then
464  if(jface.le.2) then
465  nopes=6
466  else
467  nopes=8
468  endif
469  nface=5
470  nope=15
471  elseif(lakon(nelem)(4:4).eq.'6') then
472  if(jface.le.2) then
473  nopes=3
474  else
475  nopes=4
476  endif
477  nface=5
478  nope=6
479  else
480  cycle
481  endif
482 !
483 ! determining the nodes of the face
484 !
485  if(nface.eq.4) then
486  do k=1,nopes
487  nodef(k)=kon(indexe+ifacet(k,jface))
488  enddo
489  elseif(nface.eq.5) then
490  if(nope.eq.6) then
491  do k=1,nopes
492  nodef(k)=kon(indexe+ifacew1(k,jface))
493  enddo
494  elseif(nope.eq.15) then
495  do k=1,nopes
496  nodef(k)=kon(indexe+ifacew2(k,jface))
497  enddo
498  endif
499  elseif(nface.eq.6) then
500  do k=1,nopes
501  nodef(k)=kon(indexe+ifaceq(k,jface))
502  enddo
503  endif
504 !
505 ! orthogonal projection on the element face
506 !
507  do k=1,nopes
508  do l=1,3
509  pl(l,k)=co(l,nodef(k))+vold(l,nodef(k))
510  enddo
511  enddo
512  do l=1,3
513  pproj(l)=p(l)
514  enddo
515 !
516  call attach(pl,pproj,nopes,ratio,dist,xi,et)
517 !
518  do l=1,3
519  al(l)=p(l)-pproj(l)
520  enddo
521 !
522 ! determining the jacobian vector on the surface
523 !
524  if(nopes.eq.9) then
525  call shape9q(xi,et,pl,xm,xs2,shp2,iflag)
526  elseif(nopes.eq.8) then
527  call shape8q(xi,et,pl,xm,xs2,shp2,iflag)
528  elseif(nopes.eq.4) then
529  call shape4q(xi,et,pl,xm,xs2,shp2,iflag)
530  elseif(nopes.eq.6) then
531  call shape6tri(xi,et,pl,xm,xs2,shp2,iflag)
532  elseif(nopes.eq.7) then
533  call shape7tri(xi,et,pl,xm,xs2,shp2,iflag)
534  else
535  call shape3tri(xi,et,pl,xm,xs2,shp2,iflag)
536  endif
537 !
538 ! normal on the surface
539 !
540  dm=dsqrt(xm(1)*xm(1)+xm(2)*xm(2)+xm(3)*xm(3))
541  do l=1,3
542  xn(l)=xm(l)/dm
543  enddo
544 !
545 ! distance from surface along normal (= clearance)
546 !
547  clear=al(1)*xn(1)+al(2)*xn(2)+al(3)*xn(3)
548  if((istep.eq.1).and.(iit.lt.0.d0)) then
549  if(clear.lt.0.d0) then
550  springarea(2,j)=clear
551  endif
552  endif
553  if(nmethod.eq.1) then
554  clear=clear-springarea(2,j)*(1.d0-reltime)
555  endif
556 !
557 ! check for an adjust parameter (only at the start
558 ! of the first step)
559 !
560  if((istep.eq.1).and.(iit.lt.0)) then
561  if(iset.ne.0) then
562 !
563 ! check whether node belongs to the adjust node
564 ! set
565 !
566  call nident(ialset(istartset(iset)),node,
567  & iendset(iset)-istartset(iset)+1,id)
568  if(id.gt.0) then
569  if(ialset(istartset(iset)+id-1).eq.node) then
570  do k=1,3
571  co(k,node)=co(k,node)-
572  & clear*straight(12+k,itri)
573  enddo
574  clear=0.d0
575  endif
576  endif
577  elseif(dabs(tietol(1,i)).ge.2.d0) then
578 !
579 ! adjust parameter
580 !
581  adjust=dabs(tietol(1,i))-2.d0
582  if(clear.le.adjust) then
583  do k=1,3
584  co(k,node)=co(k,node)-
585  & clear*straight(12+k,itri)
586  enddo
587  clear=0.d0
588  endif
589  endif
590  endif
591 !
592  if(int(elcon(3,1,imat)).eq.1) then
593 !
594 ! exponential overclosure
595 !
596  beta=elcon(1,1,imat)
597  c0=dlog(100.d0)/beta
598  else
599 !
600 ! linear or tabular overclosure
601 !
602  if(dabs(area).gt.0.d0) then
603  c0=c0coef*dsqrt(area)
604  else
605  c0=1.d-10
606  endif
607  endif
608  if(clear.gt.c0) then
609  isol=0
610 !
611  endif
612  endif
613 !
614  if(isol.ne.0) then
615 !
616 ! plane spring
617 !
618  ne=ne+1
619  ipkon(ne)=ifree
620  lakon(ne)='ESPRNGC '
621  ielmat(1,ne)=imat
622 !
623 ! nasym indicates whether at least one contact
624 ! spring elements exhibits friction in the present
625 ! step. If so, nasym=1, else nasym=0; nasym=1
626 ! triggers the asymmetric equation solver
627 !
628  if(ncmat_.ge.7) then
629  if(elcon(6,1,imat).gt.0) then
630  nasym=1
631  endif
632  endif
633  if(ncmat_.ge.8) then
634  if(elcon(8,1,imat).gt.0) then
635  nasym=1
636  endif
637  endif
638 !
639 c nelem=int(koncont(4,itri)/10.d0)
640 c jface=koncont(4,itri)-10*nelem
641 !
642 ! storing the area corresponding to the slave node
643 ! and the clearance if penetration takes place,
644 ! i.e. clear <0 at the start of the first step
645 !
646  springarea(1,j)=area
647 c if((istep.eq.1).and.(iit.lt.0.d0)) then
648 c if(clear.lt.0.d0) then
649 c springarea(2,j)=clear
650 c else
651 c springarea(2,j)=0.d0
652 c endif
653 c endif
654 !
655  do k=1,nopes
656  kon(ifree+k)=nodef(k)
657  enddo
658  ifree=ifree+nopes+1
659  kon(ifree)=node
660  ifree=ifree+1
661  kon(ifree)=j
662 !
663  write(lakon(ne)(8:8),'(i1)') nopes
664 !
665  indexel=indexel+1
666 !
667  if((nopes.eq.3).or.(nopes.eq.6)) then
668  if(filab(1)(3:3).eq.'C') then
669  write(27,100) setname(1:lenset)
670  100 format('*ELEMENT,TYPE=C3D4,ELSET=',a)
671  write(27,*) ne0+indexel,',',nodef(1),',',
672  & nodef(2),',',nodef(3),',',node
673  endif
674  else
675  if(filab(1)(3:3).eq.'C') then
676  write(27,101) setname(1:lenset)
677  101 format('*ELEMENT,TYPE=C3D6,ELSET=',a)
678  write(27,*) ne0+indexel,',',nodef(2),',',node,
679  & ',',
680  & nodef(3),',',nodef(1),',',node,',',nodef(4)
681  endif
682  endif
683  endif
684 !
685  enddo
686  enddo
687 !
688 ! closing the file containing the contact elements
689 !
690  if(filab(1)(3:3).eq.'C') then
691  close(27)
692  endif
693 !
694  return
subroutine near3d(xo, yo, zo, x, y, z, nx, ny, nz, xp, yp, zp, n, neighbor, k)
Definition: near3d.f:20
subroutine shape9q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape9q.f:20
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 shape7tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape7tri.f:20
subroutine stop()
Definition: stop.f:20
subroutine isortii(ix, iy, n, kflag)
Definition: isortii.f:6
static double * dist
Definition: radflowload.c:42
subroutine shape4q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape4q.f:20
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine dsort(dx, iy, n, kflag)
Definition: dsort.f:6
subroutine shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
subroutine attach(pneigh, pnode, nterms, ratio, dist, xil, etl)
Definition: attach.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)