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

Go to the source code of this file.

Functions/Subroutines

subroutine islavactive (tieset, ntie, itietri, cg, straight, co, vold, xo, yo, zo, x, y, z, nx, ny, nz, mi, imastop, nslavnode, islavnode, islavact)
 

Function/Subroutine Documentation

◆ islavactive()

subroutine islavactive ( 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, dimension(*)  mi,
integer, dimension(3,*)  imastop,
integer, dimension(*)  nslavnode,
integer, dimension(*)  islavnode,
integer, dimension(*)  islavact 
)
26 !
27  implicit none
28 !
29  character*81 tieset(3,*)
30 !
31  integer ntie,itietri(2,ntie),node,neigh(1),iflag,kneigh,i,j,k,l,
32  & isol,itri,ll,kflag,n,nx(*),ny(*),mi(*),nz(*),nstart,
33  & ifacew1(4,5),iloc,ifacew2(8,5),imastop(3,*),
34  & itriangle(100),ntriangle,ntriangle_,itriold,itrinew,id,
35  & nslavnode(*),islavnode(*),islavact(*),ifaceq(8,6),
36  & ifacet(6,4)
37 !
38  real*8 cg(3,*),straight(16,*),co(3,*),vold(0:mi(2),*),p(3),
39  & dist,xo(*),yo(*),zo(*),x(*),y(*),z(*)
40 !
41 ! nodes per face for hex elements
42 !
43  data ifaceq /4,3,2,1,11,10,9,12,
44  & 5,6,7,8,13,14,15,16,
45  & 1,2,6,5,9,18,13,17,
46  & 2,3,7,6,10,19,14,18,
47  & 3,4,8,7,11,20,15,19,
48  & 4,1,5,8,12,17,16,20/
49 !
50 ! nodes per face for tet elements
51 !
52  data ifacet /1,3,2,7,6,5,
53  & 1,2,4,5,9,8,
54  & 2,3,4,6,10,9,
55  & 1,4,3,8,10,7/
56 !
57 ! nodes per face for linear wedge elements
58 !
59  data ifacew1 /1,3,2,0,
60  & 4,5,6,0,
61  & 1,2,5,4,
62  & 2,3,6,5,
63  & 3,1,4,6/
64 !
65 ! nodes per face for quadratic wedge elements
66 !
67  data ifacew2 /1,3,2,9,8,7,0,0,
68  & 4,5,6,10,11,12,0,0,
69  & 1,2,5,4,7,14,10,13,
70  & 2,3,6,5,8,15,11,14,
71  & 3,1,4,6,9,13,12,15/
72 !
73  data iflag /2/
74 !
75 ! ***ISLAVACT***
76 !
77  iloc = 0
78 !
79  do i=1,ntie
80  if(tieset(1,i)(81:81).ne.'C') cycle
81  kneigh=1
82 !
83 ! search a master face for each slave node
84 !
85  nstart=itietri(1,i)-1
86  n=itietri(2,i)-nstart
87  if(n.lt.kneigh) kneigh=n
88  do j=1,n
89  xo(j)=cg(1,nstart+j)
90  x(j)=xo(j)
91  nx(j)=j
92  yo(j)=cg(2,nstart+j)
93  y(j)=yo(j)
94  ny(j)=j
95  zo(j)=cg(3,nstart+j)
96  z(j)=zo(j)
97  nz(j)=j
98  enddo
99  kflag=2
100  call dsort(x,nx,n,kflag)
101  call dsort(y,ny,n,kflag)
102  call dsort(z,nz,n,kflag)
103 !
104  do j=nslavnode(i)+1,nslavnode(i+1)
105  node=islavnode(j)
106 !
107  do k=1,3
108  p(k)=co(k,node)+vold(k,node)
109  enddo
110 !
111 ! determining the kneigh neighboring master contact
112 ! triangle centers of gravity
113 !
114  call near3d(xo,yo,zo,x,y,z,nx,ny,nz,p(1),p(2),p(3),
115  & n,neigh,kneigh)
116 !
117  isol=0
118 !
119  itriold=0
120  itri=neigh(1)+itietri(1,i)-1
121  ntriangle=0
122  ntriangle_=100
123 !
124  loop1: do
125  do l=1,3
126  ll=4*l-3
127  dist=straight(ll,itri)*p(1)+
128  & straight(ll+1,itri)*p(2)+
129  & straight(ll+2,itri)*p(3)+
130  & straight(ll+3,itri)
131  if(dist.gt.1.d-6) then
132  itrinew=imastop(l,itri)
133  if(itrinew.eq.0) then
134 c write(*,*) '**border reached'
135  exit loop1
136  elseif(itrinew.eq.itriold) then
137 c write(*,*) '**solution in between triangles'
138  isol=itri
139  exit loop1
140  else
141  call nident(itriangle,itrinew,ntriangle,id)
142  if(id.gt.0) then
143  if(itriangle(id).eq.itrinew) then
144 c write(*,*) '**circular path; no solution'
145  exit loop1
146  endif
147  endif
148  ntriangle=ntriangle+1
149  if(ntriangle.gt.ntriangle_) then
150 c write(*,*) '**too many iterations'
151  exit loop1
152  endif
153  do k=ntriangle,id+2,-1
154  itriangle(k)=itriangle(k-1)
155  enddo
156  itriangle(id+1)=itrinew
157  itriold=itri
158  itri=itrinew
159  cycle loop1
160  endif
161  elseif(l.eq.3) then
162 c write(*,*) '**regular solution'
163  isol=itri
164  exit loop1
165  endif
166  enddo
167  enddo loop1
168 !
169  if(isol.ne.0) then
170 ! Active node
171  islavact(j)=1
172  else
173  islavact(j)=0
174  endif
175  enddo
176  enddo
177 !
178  return
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
Hosted by OpenAircraft.com, (Michigan UAV, LLC)