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

Go to the source code of this file.

Functions/Subroutines

subroutine gentiedmpc (tieset, ntie, itietri, ipkon, kon, lakon, set, istartset, iendset, ialset, cg, straight, koncont, co, xo, yo, zo, x, y, z, nx, ny, nz, nset, ifaceslave, istartfield, iendfield, ifield, ipompc, nodempc, coefmpc, nmpc, nmpctied, mpcfree, ikmpc, ilmpc, labmpc, ithermal, tietol, cfd, ncont, imastop, ikboun, nboun, kind)
 

Function/Subroutine Documentation

◆ gentiedmpc()

subroutine gentiedmpc ( character*81, dimension(3,*)  tieset,
integer  ntie,
integer, dimension(2,ntie)  itietri,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
real*8, dimension(3,*)  cg,
real*8, dimension(16,*)  straight,
integer, dimension(4,*)  koncont,
real*8, dimension(3,*)  co,
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  nset,
integer, dimension(*)  ifaceslave,
integer, dimension(*)  istartfield,
integer, dimension(*)  iendfield,
integer, dimension(*)  ifield,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer  nmpctied,
integer  mpcfree,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
character*20, dimension(*)  labmpc,
integer, dimension(2)  ithermal,
real*8, dimension(3,*)  tietol,
integer  cfd,
integer  ncont,
integer, dimension(3,*)  imastop,
integer, dimension(*)  ikboun,
integer  nboun,
character*1  kind 
)
25 !
26 ! generates MPC's for the slave tied contact nodes
27 !
28  implicit none
29 !
30  character*1 kind
31  character*8 lakon(*)
32  character*20 labmpc(*)
33  character*81 tieset(3,*),slavset,set(*)
34 !
35  integer ntie,nset,istartset(*),iendset(*),ialset(*),
36  & itietri(2,ntie),ipkon(*),kon(*),koncont(4,*),node,
37  & neigh(1),iflag,kneigh,i,j,k,l,isol,itri,ll,kflag,n,nx(*),
38  & ny(*),nz(*),nstart,ifaceq(8,6),ifacet(6,4),nboun,
39  & ifacew1(4,5),ifacew2(8,5),nelem,jface,indexe,imastop(3,*),
40  & nnodelem,nface,nope,nodef(8),idof,kstart,kend,jstart,id,
41  & jend,ifield(*),istartfield(*),iendfield(*),ifaceslave(*),
42  & ipompc(*),nodempc(3,*),nmpc,nmpctied,mpcfree,ikmpc(*),
43  & ilmpc(*),ithermal(2),cfd,ncont,mpcfreeold,m,id1,ikboun(*),
44  & itriold,itrinew,ntriangle,ntriangle_,itriangle(100)
45 !
46  real*8 cg(3,*),straight(16,*),co(3,*),p(3),
47  & dist,xo(*),yo(*),zo(*),x(*),y(*),z(*),pl(3,9),
48  & ratio(9),xi,et,coefmpc(*),tietol(3,*),tolloc
49 !
50 ! nodes per face for hex elements
51 !
52  data ifaceq /4,3,2,1,11,10,9,12,
53  & 5,6,7,8,13,14,15,16,
54  & 1,2,6,5,9,18,13,17,
55  & 2,3,7,6,10,19,14,18,
56  & 3,4,8,7,11,20,15,19,
57  & 4,1,5,8,12,17,16,20/
58 !
59 ! nodes per face for tet elements
60 !
61  data ifacet /1,3,2,7,6,5,
62  & 1,2,4,5,9,8,
63  & 2,3,4,6,10,9,
64  & 1,4,3,8,10,7/
65 !
66 ! nodes per face for linear wedge elements
67 !
68  data ifacew1 /1,3,2,0,
69  & 4,5,6,0,
70  & 1,2,5,4,
71  & 2,3,6,5,
72  & 3,1,4,6/
73 !
74 ! nodes per face for quadratic wedge elements
75 !
76  data ifacew2 /1,3,2,9,8,7,0,0,
77  & 4,5,6,10,11,12,0,0,
78  & 1,2,5,4,7,14,10,13,
79  & 2,3,6,5,8,15,11,14,
80  & 3,1,4,6,9,13,12,15/
81 !
82 ! opening a file to store the nodes which are not connected
83 !
84  open(40,file='WarnNodeMissTiedContact.nam',status='unknown')
85  write(40,*) '*NSET,NSET=WarnNodeMissTiedContact'
86  write(*,*) '*INFO in gentiedmpc:'
87  write(*,*) ' failed nodes (if any) are stored in file'
88  write(*,*) ' WarnNodeMissTiedContact.nam'
89  write(*,*) ' This file can be loaded into'
90  write(*,*) ' an active cgx-session by typing'
91  write(*,*)
92  & ' read WarnNodeMissTiedContact.nam inp'
93  write(*,*)
94 !
95  nmpctied=nmpc
96 !
97 ! calculating a typical element size
98 !
99  tolloc=0.d0
100  do i=1,ncont
101  tolloc=tolloc+dabs(straight(1,i)*cg(1,i)+
102  & straight(2,i)*cg(2,i)+
103  & straight(3,i)*cg(3,i)+
104  & straight(4,i))
105  enddo
106  tolloc=0.025*tolloc/ncont
107 !
108 ! maximum number of neighboring master triangles for a slave node
109 !
110  kflag=2
111 !
112  do i=1,ntie
113  if(tieset(1,i)(81:81).ne.kind) cycle
114 !
115 ! determining for which dofs MPC's have to be generated
116 !
117  if(kind.eq.'T') then
118 !
119 ! thermomechanical ties or CFD
120 !
121  if(cfd.eq.1) then
122  if(ithermal(2).le.1) then
123  kstart=1
124  kend=4
125  else
126  kstart=0
127  kend=4
128  endif
129  else
130  if(ithermal(2).le.1) then
131  kstart=1
132  kend=3
133  elseif(ithermal(2).eq.2) then
134  kstart=0
135  kend=0
136  else
137  kstart=0
138  kend=3
139  endif
140  endif
141  elseif(kind.eq.'E') then
142 !
143 ! electromagnetic ties between the domains
144 !
145  if((tieset(1,i)(1:2).eq.'12').or.
146  & (tieset(1,i)(1:2).eq.'13').or.
147  & (tieset(1,i)(1:2).eq.'23')) then
148  kstart=1
149  kend=3
150  elseif((tieset(1,i)(1:2).eq.'21').or.
151  & (tieset(1,i)(1:2).eq.'31')) then
152  kstart=5
153  kend=5
154  endif
155  endif
156 !
157  iflag=0
158  kneigh=1
159  slavset=tieset(2,i)
160 !
161 ! default tolerance if none is specified
162 !
163  if(tietol(1,i).lt.1.d-10) tietol(1,i)=tolloc
164 !
165 ! determining the slave set
166 !
167  if(ifaceslave(i).eq.0) then
168  do j=1,nset
169  if(set(j).eq.slavset) then
170  exit
171  endif
172  enddo
173  jstart=istartset(j)
174  jend=iendset(j)
175  else
176  jstart=istartfield(i)
177  jend=iendfield(i)
178  endif
179 !
180  nstart=itietri(1,i)-1
181  n=itietri(2,i)-nstart
182  if(n.lt.kneigh) kneigh=n
183  do j=1,n
184  xo(j)=cg(1,nstart+j)
185  x(j)=xo(j)
186  nx(j)=j
187  yo(j)=cg(2,nstart+j)
188  y(j)=yo(j)
189  ny(j)=j
190  zo(j)=cg(3,nstart+j)
191  z(j)=zo(j)
192  nz(j)=j
193  enddo
194  call dsort(x,nx,n,kflag)
195  call dsort(y,ny,n,kflag)
196  call dsort(z,nz,n,kflag)
197 !
198  do j=jstart,jend
199  if(((ifaceslave(i).eq.0).and.(ialset(j).gt.0)).or.
200  & (ifaceslave(i).eq.1)) then
201 !
202  if(ifaceslave(i).eq.0) then
203  node=ialset(j)
204  else
205  node=ifield(j)
206  endif
207 !
208  do k=1,3
209  p(k)=co(k,node)
210  enddo
211 !
212 ! determining the kneigh neighboring master contact
213 ! triangle centers of gravity
214 !
215  call near3d(xo,yo,zo,x,y,z,nx,ny,nz,p(1),p(2),p(3),
216  & n,neigh,kneigh)
217 !
218  isol=0
219 !
220  isol=0
221 !
222  itriold=0
223  itri=neigh(1)+itietri(1,i)-1
224  ntriangle=0
225  ntriangle_=100
226 !
227  loop1: do
228  do l=1,3
229  ll=4*l-3
230  dist=straight(ll,itri)*p(1)+
231  & straight(ll+1,itri)*p(2)+
232  & straight(ll+2,itri)*p(3)+
233  & straight(ll+3,itri)
234 c if(dist.gt.1.d-6) then
235 c if(dist.gt.tietol(1,i)) then
236  if(dist.gt.tolloc) then
237  itrinew=imastop(l,itri)
238  if(itrinew.eq.0) then
239 c write(*,*) '**border reached'
240  isol=-1
241  exit loop1
242  elseif(itrinew.eq.itriold) then
243 c write(*,*) '**solution in between triangles'
244  isol=itri
245  exit loop1
246  else
247  call nident(itriangle,itrinew,ntriangle,id)
248  if(id.gt.0) then
249  if(itriangle(id).eq.itrinew) then
250 c write(*,*) '**circular path; no solution'
251  isol=-2
252  exit loop1
253  endif
254  endif
255  ntriangle=ntriangle+1
256  if(ntriangle.gt.ntriangle_) then
257 c write(*,*) '**too many iterations'
258  isol=-3
259  exit loop1
260  endif
261  do k=ntriangle,id+2,-1
262  itriangle(k)=itriangle(k-1)
263  enddo
264  itriangle(id+1)=itrinew
265  itriold=itri
266  itri=itrinew
267  cycle loop1
268  endif
269  elseif(l.eq.3) then
270 c write(*,*) '**regular solution'
271  isol=itri
272  exit loop1
273  endif
274  enddo
275  enddo loop1
276 !
277 ! if an opposite triangle is found: check the distance
278 ! perpendicular to the triangle
279 !
280  if(isol.gt.0) then
281  dist=dabs(straight(13,itri)*p(1)+
282  & straight(14,itri)*p(2)+
283  & straight(15,itri)*p(3)+
284  & straight(16,itri))
285  if(dist.gt.tietol(1,i)) isol=0
286  endif
287 !
288  if(isol.le.0) then
289 !
290 ! no MPC is generated
291 !
292  write(*,*) '*WARNING in gentiedmpc: no tied MPC'
293  write(*,*) ' generated for node ',node
294  if(isol.eq.0) then
295  write(*,*) ' master face too far away'
296  write(*,*) ' distance: ',dist
297  write(*,*) ' tolerance: ',tietol(1,i)
298  else
299  write(*,*) ' no corresponding master face'
300  write(*,*) ' found; tolerance: ',
301  & tolloc
302 c & tietol(1,i)
303  endif
304  write(40,*) node
305  else
306 !
307  nelem=int(koncont(4,itri)/10.d0)
308  jface=koncont(4,itri)-10*nelem
309 !
310  indexe=ipkon(nelem)
311  if(lakon(nelem)(4:4).eq.'2') then
312  nnodelem=8
313  nface=6
314  elseif(lakon(nelem)(4:4).eq.'8') then
315  nnodelem=4
316  nface=6
317  elseif(lakon(nelem)(4:5).eq.'10') then
318  nnodelem=6
319  nface=4
320  elseif(lakon(nelem)(4:4).eq.'4') then
321  nnodelem=3
322  nface=4
323  elseif(lakon(nelem)(4:5).eq.'15') then
324  if(jface.le.2) then
325  nnodelem=6
326  else
327  nnodelem=8
328  endif
329  nface=5
330  nope=15
331  elseif(lakon(nelem)(4:4).eq.'6') then
332  if(jface.le.2) then
333  nnodelem=3
334  else
335  nnodelem=4
336  endif
337  nface=5
338  nope=6
339  else
340  cycle
341  endif
342 !
343 ! determining the nodes of the face
344 !
345  if(nface.eq.4) then
346  do k=1,nnodelem
347  nodef(k)=kon(indexe+ifacet(k,jface))
348  enddo
349  elseif(nface.eq.5) then
350  if(nope.eq.6) then
351  do k=1,nnodelem
352  nodef(k)=kon(indexe+ifacew1(k,jface))
353  enddo
354  elseif(nope.eq.15) then
355  do k=1,nnodelem
356  nodef(k)=kon(indexe+ifacew2(k,jface))
357  enddo
358  endif
359  elseif(nface.eq.6) then
360  do k=1,nnodelem
361  nodef(k)=kon(indexe+ifaceq(k,jface))
362  enddo
363  endif
364 !
365 ! attaching the node with coordinates in p
366 ! to the face
367 !
368  do k=1,nnodelem
369  do l=1,3
370  pl(l,k)=co(l,nodef(k))
371  enddo
372  enddo
373  call attach(pl,p,nnodelem,ratio,dist,xi,et)
374 !
375 ! adjusting the coordinates of the node (only
376 ! if the user did not specify ADJUST=NO on the
377 ! *TIE card)
378 !
379  if(tietol(2,i).gt.0.d0) then
380  do k=1,3
381  co(k,node)=p(k)
382  enddo
383  endif
384 !
385 ! generating MPC's
386 !
387  do l=kstart,kend
388  idof=8*(node-1)+l
389  call nident(ikmpc,idof,nmpc,id)
390  if(id.gt.0) then
391  if(ikmpc(id).eq.idof) then
392  write(*,*) '*WARNING in gentiedmpc:'
393  write(*,*) ' DOF ',l,' of node ',
394  & node,' is not active;'
395  write(*,*) ' no tied constraint ',
396  & 'is generated'
397  write(40,*) node
398  cycle
399  endif
400  endif
401 !
402  call nident(ikboun,idof,nboun,id1)
403  if(id1.gt.0) then
404  if(ikboun(id1).eq.idof) then
405  write(*,*) '*WARNING in gentiedmpc:'
406  write(*,*) ' DOF ',l,' of node ',
407  & node,' is not active;'
408  write(*,*) ' no tied constraint ',
409  & 'is generated'
410  write(40,*) node
411  cycle
412  endif
413  endif
414 !
415  nmpc=nmpc+1
416  labmpc(nmpc)=' '
417  ipompc(nmpc)=mpcfree
418 !
419 ! updating ikmpc and ilmpc
420 !
421  do m=nmpc,id+2,-1
422  ikmpc(m)=ikmpc(m-1)
423  ilmpc(m)=ilmpc(m-1)
424  enddo
425  ikmpc(id+1)=idof
426  ilmpc(id+1)=nmpc
427 !
428  nodempc(1,mpcfree)=node
429  nodempc(2,mpcfree)=l
430  coefmpc(mpcfree)=1.d0
431  mpcfree=nodempc(3,mpcfree)
432  if(mpcfree.eq.0) then
433  write(*,*)
434  & '*ERROR in gentiedmpc: increase memmpc_'
435  call exit(201)
436  endif
437  do k=1,nnodelem
438  nodempc(1,mpcfree)=nodef(k)
439  nodempc(2,mpcfree)=l
440  coefmpc(mpcfree)=-ratio(k)
441  mpcfreeold=mpcfree
442  mpcfree=nodempc(3,mpcfree)
443  if(mpcfree.eq.0) then
444  write(*,*)
445  & '*ERROR in gentiedmpc: increase memmpc_'
446  call exit(201)
447  endif
448  enddo
449  nodempc(3,mpcfreeold)=0
450  enddo
451 !
452  endif
453 !
454  else
455  node=ialset(j-2)
456  do
457  node=node-ialset(j)
458  if(node.ge.ialset(j-1)) exit
459 !
460  do k=1,3
461  p(k)=co(k,node)
462  enddo
463 !
464 ! determining the kneigh neighboring master contact
465 ! triangle centers of gravity
466 !
467  call near3d(xo,yo,zo,x,y,z,nx,ny,nz,p(1),p(2),p(3),
468  & n,neigh,kneigh)
469 !
470  isol=0
471 !
472  isol=0
473 !
474  itriold=0
475  itri=neigh(1)+itietri(1,i)-1
476  ntriangle=0
477  ntriangle_=100
478 !
479  loop2: do
480  do l=1,3
481  ll=4*l-3
482  dist=straight(ll,itri)*p(1)+
483  & straight(ll+1,itri)*p(2)+
484  & straight(ll+2,itri)*p(3)+
485  & straight(ll+3,itri)
486 c if(dist.gt.1.d-6) then
487 c if(dist.gt.tietol(1,i)) then
488  if(dist.gt.tolloc) then
489  itrinew=imastop(l,itri)
490  if(itrinew.eq.0) then
491 c write(*,*) '**border reached'
492  isol=-1
493  exit loop2
494  elseif(itrinew.eq.itriold) then
495 c write(*,*) '**solution in between triangles'
496  isol=itri
497  exit loop2
498  else
499  call nident(itriangle,itrinew,ntriangle,id)
500  if(id.gt.0) then
501  if(itriangle(id).eq.itrinew) then
502 c write(*,*) '**circular path; no solution'
503  isol=-2
504  exit loop2
505  endif
506  endif
507  ntriangle=ntriangle+1
508  if(ntriangle.gt.ntriangle_) then
509 c write(*,*) '**too many iterations'
510  isol=-3
511  exit loop2
512  endif
513  do k=ntriangle,id+2,-1
514  itriangle(k)=itriangle(k-1)
515  enddo
516  itriangle(id+1)=itrinew
517  itriold=itri
518  itri=itrinew
519  cycle loop2
520  endif
521  elseif(l.eq.3) then
522 c write(*,*) '**regular solution'
523  isol=itri
524  exit loop2
525  endif
526  enddo
527  enddo loop2
528 !
529 ! if an opposite triangle is found: check the distance
530 ! perpendicular to the triangle
531 !
532  if(isol.gt.0) then
533  dist=dabs(straight(13,itri)*p(1)+
534  & straight(14,itri)*p(2)+
535  & straight(15,itri)*p(3)+
536  & straight(16,itri))
537  if(dist.gt.tietol(1,i)) isol=0
538  endif
539 !
540 ! check whether distance is larger than tietol(1,i):
541 ! no element is generated
542 !
543  if(isol.eq.0) then
544 !
545 ! no MPC is generated
546 !
547  write(*,*) '*WARNING in gentiedmpc: no tied MPC'
548  write(*,*) ' generated for node ',node
549  if(isol.eq.0) then
550  write(*,*) ' master face too far away'
551  write(*,*) ' distance: ',dist
552  write(*,*) ' tolerance: ',tietol(1,i)
553  else
554  write(*,*) ' no corresponding master face'
555  write(*,*) ' found; tolerance: ',
556  & tolloc
557 c & tietol(1,i)
558  endif
559  write(40,*) node
560  else
561 !
562  nelem=int(koncont(4,itri)/10.d0)
563  jface=koncont(4,itri)-10*nelem
564 !
565  indexe=ipkon(nelem)
566  if(lakon(nelem)(4:4).eq.'2') then
567  nnodelem=8
568  nface=6
569  elseif(lakon(nelem)(4:4).eq.'8') then
570  nnodelem=4
571  nface=6
572  elseif(lakon(nelem)(4:5).eq.'10') then
573  nnodelem=6
574  nface=4
575  elseif(lakon(nelem)(4:4).eq.'4') then
576  nnodelem=3
577  nface=4
578  elseif(lakon(nelem)(4:5).eq.'15') then
579  if(jface.le.2) then
580  nnodelem=6
581  else
582  nnodelem=8
583  endif
584  nface=5
585  nope=15
586  elseif(lakon(nelem)(4:4).eq.'6') then
587  if(jface.le.2) then
588  nnodelem=3
589  else
590  nnodelem=4
591  endif
592  nface=5
593  nope=6
594  else
595  cycle
596  endif
597 !
598 ! determining the nodes of the face
599 !
600  if(nface.eq.4) then
601  do k=1,nnodelem
602  nodef(k)=kon(indexe+ifacet(k,jface))
603  enddo
604  elseif(nface.eq.5) then
605  if(nope.eq.6) then
606  do k=1,nnodelem
607  nodef(k)=kon(indexe+ifacew1(k,jface))
608  enddo
609  elseif(nope.eq.15) then
610  do k=1,nnodelem
611  nodef(k)=kon(indexe+ifacew2(k,jface))
612  enddo
613  endif
614  elseif(nface.eq.6) then
615  do k=1,nnodelem
616  nodef(k)=kon(indexe+ifaceq(k,jface))
617  enddo
618  endif
619 !
620 ! attaching the node with coordinates in p
621 ! to the face
622 !
623  do k=1,nnodelem
624  do l=1,3
625  pl(l,k)=co(l,nodef(k))
626  enddo
627  enddo
628  call attach(pl,p,nnodelem,ratio,dist,xi,et)
629 !
630 ! adjusting the coordinates of the node (only
631 ! if the user did not specify ADJUST=NO on the
632 ! *TIE card)
633 !
634  if(tietol(2,i).gt.0.d0) then
635  do k=1,3
636  co(k,node)=p(k)
637  enddo
638  endif
639 !
640 ! generating MPC's
641 !
642  do l=kstart,kend
643  idof=8*(node-1)+l
644  call nident(ikmpc,idof,nmpc,id)
645  if(id.gt.0) then
646  if(ikmpc(id).eq.idof) then
647  write(*,*) '*WARNING in gentiedmpc:'
648  write(*,*) ' DOF ',l,' of node ',
649  & node,' is not active;'
650  write(*,*) ' no tied constraint ',
651  & 'is generated'
652  write(40,*) node
653  cycle
654  endif
655  endif
656 !
657  call nident(ikboun,idof,nboun,id1)
658  if(id1.gt.0) then
659  if(ikboun(id1).eq.idof) then
660  write(*,*) '*WARNING in gentiedmpc:'
661  write(*,*) ' DOF ',l,' of node ',
662  & node,' is not active;'
663  write(*,*) ' no tied constraint ',
664  & 'is generated'
665  write(40,*) node
666  cycle
667  endif
668  endif
669 !
670  nmpc=nmpc+1
671  labmpc(nmpc)=' '
672  ipompc(nmpc)=mpcfree
673 !
674 ! updating ikmpc and ilmpc
675 !
676  do m=nmpc,id+2,-1
677  ikmpc(m)=ikmpc(m-1)
678  ilmpc(m)=ilmpc(m-1)
679  enddo
680  ikmpc(id+1)=idof
681  ilmpc(id+1)=nmpc
682 !
683  nodempc(1,mpcfree)=node
684  nodempc(2,mpcfree)=l
685  coefmpc(mpcfree)=1.d0
686  mpcfree=nodempc(3,mpcfree)
687  if(mpcfree.eq.0) then
688  write(*,*)
689  & '*ERROR in gentiedmpc: increase memmpc_'
690  call exit(201)
691  endif
692  do k=1,nnodelem
693  nodempc(1,mpcfree)=nodef(k)
694  nodempc(2,mpcfree)=l
695  coefmpc(mpcfree)=-ratio(k)
696  mpcfreeold=mpcfree
697  mpcfree=nodempc(3,mpcfree)
698  if(mpcfree.eq.0) then
699  write(*,*)
700  & '*ERROR in gentiedmpc: increase memmpc_'
701  call exit(201)
702  endif
703  enddo
704  nodempc(3,mpcfreeold)=0
705  enddo
706  endif
707 !
708  enddo
709  endif
710  enddo
711  enddo
712 !
713 ! number of tied MPC's
714 !
715  nmpctied=nmpc-nmpctied
716 !
717  close(40)
718 !
719  return
720 !
subroutine near3d(xo, yo, zo, x, y, z, nx, ny, nz, xp, yp, zp, n, neighbor, k)
Definition: near3d.f:20
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
subroutine attach(pneigh, pnode, nterms, ratio, dist, xil, etl)
Definition: attach.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)