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

Go to the source code of this file.

Functions/Subroutines

subroutine normalsoninterface (istartset, iendset, ialset, imast, ipkon, kon, lakon, imastnode, nmastnode, xmastnor, co)
 

Function/Subroutine Documentation

◆ normalsoninterface()

subroutine normalsoninterface ( integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer  imast,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
integer, dimension(*)  imastnode,
integer  nmastnode,
real*8, dimension(3,*)  xmastnor,
real*8, dimension(3,*)  co 
)
22 !
23  logical exist
24 !
25  character*8 lakon(*)
26 !
27  integer istartset(*),iendset(*),ialset(*),imast,j,ifacem,
28  & nelemm,jfacem,indexe,ipkon(*),kon(*),nopem,node,
29  & ifaceq(8,6),ifacet(6,4),ifacew1(4,5),ifacew2(8,5),
30  & nmastnode,id,indexnode(9),konl(26),imastnode(*)
31 !
32  real*8 xmastnor(3,*),xsj2(3),shp2(7,9),xs2(3,2),xi,et,dd,
33  & xquad(2,9),xtri(2,7),xl2m(3,9),co(3,*)
34 !
35 ! nodes per face for hex elements
36 !
37  data ifaceq /4,3,2,1,11,10,9,12,
38  & 5,6,7,8,13,14,15,16,
39  & 1,2,6,5,9,18,13,17,
40  & 2,3,7,6,10,19,14,18,
41  & 3,4,8,7,11,20,15,19,
42  & 4,1,5,8,12,17,16,20/
43 !
44 ! nodes per face for tet elements
45 !
46  data ifacet /1,3,2,7,6,5,
47  & 1,2,4,5,9,8,
48  & 2,3,4,6,10,9,
49  & 1,4,3,8,10,7/
50 !
51 ! nodes per face for linear wedge elements
52 !
53  data ifacew1 /1,3,2,0,
54  & 4,5,6,0,
55  & 1,2,5,4,
56  & 2,3,6,5,
57  & 3,1,4,6/
58 !
59 ! nodes per face for quadratic wedge elements
60 !
61  data ifacew2 /1,3,2,9,8,7,0,0,
62  & 4,5,6,10,11,12,0,0,
63  & 1,2,5,4,7,14,10,13,
64  & 2,3,6,5,8,15,11,14,
65  & 3,1,4,6,9,13,12,15/
66 !
67 ! new added data for the local coodinates for nodes
68 !
69  data xquad /-1.d0,-1.d0,
70  & 1.d0,-1.d0,
71  & 1.d0,1.d0,
72  & -1.d0,1.d0,
73  & 0.d0,-1.d0,
74  & 1.d0,0.d0,
75  & 0.d0,1.d0,
76  & -1.d0,0.d0,
77  & 0.d0,0.d0/
78 !
79  data xtri /0.d0,0.d0,
80  & 1.d0,0.d0,
81  & 0.d0,1.d0,
82  & .5d0,0.d0,
83  & .5d0,.5d0,
84  & 0.d0,.5d0,
85  & 0.333333333333333d0,0.333333333333333d0/
86 !
87  data iflag /2/
88 !
89 ! creates the normals for the interface condition A.n=0
90 ! for electromagnetic calculations
91 !
92  nmastnode=0
93 !
94  do j=istartset(imast),iendset(imast)
95 !
96 ! create imastnode, and nmastnode (number of nodes in the
97 ! master surface
98 !
99  ifacem=ialset(j)
100  nelemm=int(ifacem/10)
101  jfacem=ifacem - nelemm*10
102  indexe=ipkon(nelemm)
103 !
104  if(lakon(nelemm)(4:5).eq.'20') then
105  nopem=8
106  elseif(lakon(nelemm)(4:4).eq.'8') then
107  nopem=4
108  elseif(lakon(nelemm)(4:5).eq.'10') then
109  nopem=6
110  elseif(lakon(nelemm)(4:4).eq.'4') then
111  nopem=3
112  endif
113 !
114  if(lakon(nelemm)(4:4).eq.'6') then
115  if(jfacem.le.2) then
116  nopem=3
117  else
118  nopem=4
119  endif
120  endif
121  if(lakon(nelemm)(4:5).eq.'15') then
122  if(jfacem.le.2) then
123  nopem=6
124  else
125  nopem=8
126  endif
127  endif
128 !
129  do l=1,nopem
130  if((lakon(nelemm)(4:4).eq.'2').or.
131  & (lakon(nelemm)(4:4).eq.'8')) then
132  node=kon(indexe+ifaceq(l,jfacem))
133  elseif((lakon(nelemm)(4:4).eq.'4').or.
134  & (lakon(nelemm)(4:5).eq.'10')) then
135  node=kon(indexe+ifacet(l,jfacem))
136  elseif(lakon(nelemm)(4:4).eq.'6') then
137  node=kon(indexe+ifacew1(l,jfacem))
138  elseif(lakon(nelemm)(4:5).eq.'15') then
139  node=kon(indexe+ifacew2(l,jfacem))
140  endif
141  call nident(imastnode,node,
142  & nmastnode,id)
143  exist=.false.
144  if(id.gt.0) then
145  if(imastnode(id).eq.node) then
146  exist=.true.
147  endif
148  endif
149  if(exist) cycle
150  nmastnode=nmastnode+1
151  do k=nmastnode,id+2,-1
152  imastnode(k)=imastnode(k-1)
153  enddo
154  imastnode(id+1)=node
155  enddo
156 !
157  enddo
158 !
159 ! calculate the mean normal in the master nodes
160 !
161  do j=istartset(imast),iendset(imast)
162 !
163  ifacem=ialset(j)
164  nelemm=int(ifacem/10)
165  jfacem=ifacem - nelemm*10
166  indexe=ipkon(nelemm)
167 !
168 ! nopem: # of nodes in the master face
169 ! nope: # of nodes in the element
170 !
171  if(lakon(nelemm)(4:5).eq.'8R') then
172  nopem=4
173  nope=8
174  elseif(lakon(nelemm)(4:4).eq.'8') then
175  nopem=4
176  nope=8
177  elseif(lakon(nelemm)(4:5).eq.'20') then
178  nopem=8
179  nope=20
180  elseif(lakon(nelemm)(4:5).eq.'10') then
181  nopem=6
182  nope=10
183  elseif(lakon(nelemm)(4:4).eq.'4') then
184  nopem=3
185  nope=4
186 !
187 ! treatment of wedge faces
188 !
189  elseif(lakon(nelemm)(4:4).eq.'6') then
190  nope=6
191  if(jfacem.le.2) then
192  nopem=3
193  else
194  nopem=4
195  endif
196  elseif(lakon(nelemm)(4:5).eq.'15') then
197  nope=15
198  if(jfacem.le.2) then
199  nopem=6
200  else
201  nopem=8
202  endif
203  endif
204 !
205 ! actual position of the nodes belonging to the
206 ! master surface
207 !
208  do k=1,nope
209  konl(k)=kon(ipkon(nelemm)+k)
210  enddo
211 !
212  if((nope.eq.20).or.(nope.eq.8)) then
213  do m=1,nopem
214  do k=1,3
215  xl2m(k,m)=co(k,konl(ifaceq(m,jfacem)))
216  enddo
217  enddo
218  elseif((nope.eq.10).or.(nope.eq.4))
219  & then
220  do m=1,nopem
221  do k=1,3
222  xl2m(k,m)=co(k,konl(ifacet(m,jfacem)))
223  enddo
224  enddo
225  elseif(nope.eq.15) then
226  do m=1,nopem
227  do k=1,3
228  xl2m(k,m)=co(k,konl(ifacew2(m,jfacem)))
229  enddo
230  enddo
231  else
232  do m=1,nopem
233  do k=1,3
234  xl2m(k,m)=co(k,konl(ifacew1(m,jfacem)))
235  enddo
236  enddo
237  endif
238 
239 ! calculate the normal vector in the nodes belonging to the master surface
240 !
241  if(nopem.eq.8) then
242  do m=1,nopem
243  xi=xquad(1,m)
244  et=xquad(2,m)
245  call shape8q(xi,et,xl2m,xsj2,xs2,shp2,iflag)
246  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
247  & + xsj2(3)*xsj2(3))
248  xsj2(1)=xsj2(1)/dd
249  xsj2(2)=xsj2(2)/dd
250  xsj2(3)=xsj2(3)/dd
251 !
252  if(nope.eq.20) then
253  node=konl(ifaceq(m,jfacem))
254  elseif(nope.eq.15) then
255  node=konl(ifacew2(m,jfacem))
256  endif
257 !
258  call nident(imastnode,node,
259  & nmastnode,id)
260  index1=id
261  indexnode(m)=index1
262  xmastnor(1,index1)=xmastnor(1,index1)
263  & +xsj2(1)
264  xmastnor(2,index1)=xmastnor(2,index1)
265  & +xsj2(2)
266  xmastnor(3,index1)=xmastnor(3,index1)
267  & +xsj2(3)
268  enddo
269  elseif(nopem.eq.4) then
270  do m=1,nopem
271  xi=xquad(1,m)
272  et=xquad(2,m)
273  call shape4q(xi,et,xl2m,xsj2,xs2,shp2,iflag)
274  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
275  & + xsj2(3)*xsj2(3))
276  xsj2(1)=xsj2(1)/dd
277  xsj2(2)=xsj2(2)/dd
278  xsj2(3)=xsj2(3)/dd
279 !
280  if(nope.eq.8) then
281  node=konl(ifaceq(m,jfacem))
282  elseif(nope.eq.6) then
283  node=konl(ifacew1(m,jfacem))
284  endif
285 !
286  call nident(imastnode,node,
287  & nmastnode,id)
288 !
289  index1=id
290  indexnode(m)=index1
291  xmastnor(1,index1)=xmastnor(1,index1)
292  & +xsj2(1)
293  xmastnor(2,index1)=xmastnor(2,index1)
294  & +xsj2(2)
295  xmastnor(3,index1)=xmastnor(3,index1)
296  & +xsj2(3)
297  enddo
298  elseif(nopem.eq.6) then
299  do m=1,nopem
300  xi=xtri(1,m)
301  et=xtri(2,m)
302  call shape6tri(xi,et,xl2m,xsj2,xs2,shp2,iflag)
303  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
304  & + xsj2(3)*xsj2(3))
305  xsj2(1)=xsj2(1)/dd
306  xsj2(2)=xsj2(2)/dd
307  xsj2(3)=xsj2(3)/dd
308 !
309  if(nope.eq.10) then
310  node=konl(ifacet(m,jfacem))
311  elseif(nope.eq.15) then
312  node=konl(ifacew2(m,jfacem))
313  endif
314 !
315  call nident(imastnode,node,
316  & nmastnode,id)
317  index1=id
318  indexnode(m)=index1
319  xmastnor(1,index1)=xmastnor(1,index1)
320  & +xsj2(1)
321  xmastnor(2,index1)=xmastnor(2,index1)
322  & +xsj2(2)
323  xmastnor(3,index1)=xmastnor(3,index1)
324  & +xsj2(3)
325  enddo
326  else
327  do m=1,nopem
328  xi=xtri(1,m)
329  et=xtri(2,m)
330  call shape3tri(xi,et,xl2m,xsj2,xs2,shp2,iflag)
331  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
332  & + xsj2(3)*xsj2(3))
333  xsj2(1)=xsj2(1)/dd
334  xsj2(2)=xsj2(2)/dd
335  xsj2(3)=xsj2(3)/dd
336 !
337  if(nope.eq.6) then
338  node=konl(ifacew1(m,jfacem))
339  elseif(nope.eq.4) then
340  node=konl(ifacet(m,jfacem))
341  endif
342 !
343  call nident(imastnode,node,
344  & nmastnode,id)
345  index1=id
346  indexnode(m)=index1
347  xmastnor(1,id)=xmastnor(1,index1)
348  & +xsj2(1)
349  xmastnor(2,id)=xmastnor(2,index1)
350  & +xsj2(2)
351  xmastnor(3,id)=xmastnor(3,index1)
352  & +xsj2(3)
353  enddo
354  endif
355  enddo
356 !
357 ! normalizing the normals
358 !
359  do l=1,nmastnode
360  dd=dsqrt(xmastnor(1,l)**2+xmastnor(2,l)**2+
361  & xmastnor(3,l)**2)
362  do m=1,3
363  xmastnor(m,l)=xmastnor(m,l)/dd
364  enddo
365  enddo
366 !
367  return
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 shape4q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape4q.f:20
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)