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

Go to the source code of this file.

Functions/Subroutines

subroutine findsurface_se (nodface, ipoface, ne, ipkon, lakon, kon, konfa, ipkonfa, nk, lakonfa, nsurfs)
 

Function/Subroutine Documentation

◆ findsurface_se()

subroutine findsurface_se ( integer, dimension(5,*)  nodface,
integer, dimension(*)  ipoface,
integer  ne,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer, dimension(*)  kon,
integer, dimension(*)  konfa,
integer, dimension(*)  ipkonfa,
integer  nk,
character*8, dimension(*)  lakonfa,
integer  nsurfs 
)
21 !
22  implicit none
23 !
24 ! finds the external surface of the structure
25 !
26  character*8 lakon(*),lakonfa(*)
27 !
28  integer ipkon(*),kon(*),ne,nodface(5,*),ipoface(*),nk,nopem,m,
29  & ithree,ifour,ifaceq(8,6),konfa(*),ipkonfa(*),nelemm,jfacem,
30  & ifacet(6,4),ifacew2(8,5),ifree,ifreenew,index,indexold,kflag,
31  & i,j,k,nodes(4),iaux,indexe,konl(26),nope,nsurfs,
32  & ifacew1(4,5)
33 !
34 ! nodes belonging to the element faces
35 !
36  data ifaceq /4,3,2,1,11,10,9,12,
37  & 5,6,7,8,13,14,15,16,
38  & 1,2,6,5,9,18,13,17,
39  & 2,3,7,6,10,19,14,18,
40  & 3,4,8,7,11,20,15,19,
41  & 4,1,5,8,12,17,16,20/
42  data ifacet /1,3,2,7,6,5,
43  & 1,2,4,5,9,8,
44  & 2,3,4,6,10,9,
45  & 1,4,3,8,10,7/
46  data ifacew1 /1,3,2,0,
47  & 4,5,6,0,
48  & 1,2,5,4,
49  & 2,3,6,5,
50  & 3,1,4,6/
51  data ifacew2 /1,3,2,9,8,7,0,0,
52  & 4,5,6,10,11,12,0,0,
53  & 1,2,5,4,7,14,10,13,
54  & 2,3,6,5,8,15,11,14,
55  & 3,1,4,6,9,13,12,15/
56 !
57  kflag=1
58  ithree=3
59  ifour=4
60 !
61 ! determining the external element faces of the electromagnetic mesh
62 ! the faces are catalogued by the three lowes nodes numbers
63 ! in ascending order. ipoface(i) points to a face for which
64 ! node i is the lowest node and nodface(1,ipoface(i)) and
65 ! nodface(2,ipoface(i)) are the next lower ones.
66 ! nodface(3,ipoface(i)) contains the element number,
67 ! nodface(4,ipoface(i)) the face number and nodface(5,ipoface(i))
68 ! is a pointer to the next surface for which node i is the
69 ! lowest node; if there are no more such surfaces the pointer
70 ! has the value zero
71 ! An external element face is one which belongs to one element
72 ! only
73 !
74  ifree=1
75  do i=1,6*ne-1
76  nodface(5,i)=i+1
77  enddo
78  do i=1,ne
79  if(ipkon(i).lt.0) cycle
80  if(lakon(i)(1:3).ne.'C3D') cycle
81  indexe=ipkon(i)
82  if((lakon(i)(4:4).eq.'2').or.(lakon(i)(4:4).eq.'8')) then
83  do j=1,6
84  do k=1,4
85  nodes(k)=kon(indexe+ifaceq(k,j))
86  enddo
87  call isortii(nodes,iaux,ifour,kflag)
88  indexold=0
89  index=ipoface(nodes(1))
90  do
91 !
92 ! adding a surface which has not been
93 ! catalogued so far
94 !
95  if(index.eq.0) then
96  ifreenew=nodface(5,ifree)
97  nodface(1,ifree)=nodes(2)
98  nodface(2,ifree)=nodes(3)
99  nodface(3,ifree)=i
100  nodface(4,ifree)=j
101  nodface(5,ifree)=ipoface(nodes(1))
102  ipoface(nodes(1))=ifree
103  ifree=ifreenew
104  exit
105  endif
106 !
107 ! removing a surface which has already
108 ! been catalogued
109 !
110  if((nodface(1,index).eq.nodes(2)).and.
111  & (nodface(2,index).eq.nodes(3))) then
112  if(indexold.eq.0) then
113  ipoface(nodes(1))=nodface(5,index)
114  else
115  nodface(5,indexold)=nodface(5,index)
116  endif
117  nodface(5,index)=ifree
118  ifree=index
119  exit
120  endif
121  indexold=index
122  index=nodface(5,index)
123  enddo
124  enddo
125  elseif((lakon(i)(4:4).eq.'4').or.(lakon(i)(4:5).eq.'10')) then
126  do j=1,4
127  do k=1,3
128  nodes(k)=kon(indexe+ifacet(k,j))
129  enddo
130  call isortii(nodes,iaux,ithree,kflag)
131  indexold=0
132  index=ipoface(nodes(1))
133  do
134 !
135 ! adding a surface which has not been
136 ! catalogued so far
137 !
138  if(index.eq.0) then
139  ifreenew=nodface(5,ifree)
140  nodface(1,ifree)=nodes(2)
141  nodface(2,ifree)=nodes(3)
142  nodface(3,ifree)=i
143  nodface(4,ifree)=j
144  nodface(5,ifree)=ipoface(nodes(1))
145  ipoface(nodes(1))=ifree
146  ifree=ifreenew
147  exit
148  endif
149 !
150 ! removing a surface which has already
151 ! been catalogued
152 !
153  if((nodface(1,index).eq.nodes(2)).and.
154  & (nodface(2,index).eq.nodes(3))) then
155  if(indexold.eq.0) then
156  ipoface(nodes(1))=nodface(5,index)
157  else
158  nodface(5,indexold)=nodface(5,index)
159  endif
160  nodface(5,index)=ifree
161  ifree=index
162  exit
163  endif
164  indexold=index
165  index=nodface(5,index)
166  enddo
167  enddo
168  else
169  do j=1,5
170  if(j.le.2) then
171  do k=1,3
172  nodes(k)=kon(indexe+ifacew2(k,j))
173  enddo
174  call isortii(nodes,iaux,ithree,kflag)
175  else
176  do k=1,4
177  nodes(k)=kon(indexe+ifacew2(k,j))
178  enddo
179  call isortii(nodes,iaux,ifour,kflag)
180  endif
181  indexold=0
182  index=ipoface(nodes(1))
183  do
184 !
185 ! adding a surface which has not been
186 ! catalogued so far
187 !
188  if(index.eq.0) then
189  ifreenew=nodface(5,ifree)
190  nodface(1,ifree)=nodes(2)
191  nodface(2,ifree)=nodes(3)
192  nodface(3,ifree)=i
193  nodface(4,ifree)=j
194  nodface(5,ifree)=ipoface(nodes(1))
195  ipoface(nodes(1))=ifree
196  ifree=ifreenew
197  exit
198  endif
199 !
200 ! removing a surface which has already
201 ! been catalogued
202 !
203  if((nodface(1,index).eq.nodes(2)).and.
204  & (nodface(2,index).eq.nodes(3))) then
205  if(indexold.eq.0) then
206  ipoface(nodes(1))=nodface(5,index)
207  else
208  nodface(5,indexold)=nodface(5,index)
209  endif
210  nodface(5,index)=ifree
211  ifree=index
212  exit
213  endif
214  indexold=index
215  index=nodface(5,index)
216  enddo
217  enddo
218  endif
219  enddo
220 !
221 ! Create fields konfa and ipkonfa for calculation of normals of shell
222 ! elements
223 !
224  nsurfs=0
225  ifree=0
226  do j=1,nk
227 !
228  if(ipoface(j).eq.0) cycle
229  indexe=ipoface(j)
230 !
231  do
232  nsurfs=nsurfs+1
233  nelemm=nodface(3,indexe)
234  jfacem=nodface(4,indexe)
235 !
236 ! treatment of hexahedral elements
237 !
238  if(lakon(nelemm)(4:4).eq.'8') then
239  nopem=4
240  nope=8
241  elseif(lakon(nelemm)(4:5).eq.'20') then
242  nopem=8
243  nope=20
244 !
245 ! treatment of tethrahedral elements
246 !
247  elseif(lakon(nelemm)(4:5).eq.'10') then
248  nopem=6
249  nope=10
250  elseif(lakon(nelemm)(4:4).eq.'4') then
251  nopem=3
252  nope=4
253 !
254 ! treatment of wedge faces
255 !
256  elseif(lakon(nelemm)(4:4).eq.'6') then
257  nope=6
258  if(jfacem.le.2) then
259  nopem=3
260  else
261  nopem=4
262  endif
263  elseif(lakon(nelemm)(4:5).eq.'15') then
264  nope=15
265  if(jfacem.le.2) then
266  nopem=6
267  else
268  nopem=8
269  endif
270  endif
271 !
272 ! actual position of the nodes
273 !
274  do k=1,nope
275  konl(k)=kon(ipkon(nelemm)+k)
276  enddo
277 !
278 ! quadratic quad shell element
279 !
280  if((nope.eq.20).or.
281  & ((nope.eq.15).and.(jfacem.gt.2))) then
282  lakonfa(nsurfs)='S8'
283  ipkonfa(nsurfs)=ifree
284  do m=1,nopem
285  if(nope.eq.20) then
286  konfa(ifree+m)=konl(ifaceq(m,jfacem))
287  else
288  konfa(ifree+m)=konl(ifacew2(m,jfacem))
289  endif
290  enddo
291  ifree=ifree+nopem
292 !
293 ! linear quad shell element
294 !
295  elseif((nope.eq.8).or.
296  & ((nope.eq.6).and.(jfacem.gt.2))) then
297  lakonfa(nsurfs)='S4'
298  ipkonfa(nsurfs)=ifree
299  do m=1,nopem
300  if(nope.eq.8) then
301  konfa(ifree+m)=konl(ifaceq(m,jfacem))
302  else
303  konfa(ifree+m)=konl(ifacew1(m,jfacem))
304  endif
305  enddo
306  ifree=ifree+nopem
307 !
308 ! quadratic tri shell element
309 !
310  elseif((nope.eq.10).or.
311  & ((nope.eq.15).and.(jfacem.le.2))) then
312  lakonfa(nsurfs)='S6'
313  ipkonfa(nsurfs)=ifree
314  do m=1,nopem
315  if(nope.eq.10) then
316  konfa(ifree+m)=konl(ifacet(m,jfacem))
317  else
318  konfa(ifree+m)=konl(ifacew2(m,jfacem))
319  endif
320  enddo
321  ifree=ifree+nopem
322 !
323 ! linear tri shell element
324 !
325  elseif((nope.eq.4).or.
326  & ((nope.eq.6).and.(jfacem.le.2))) then
327  lakonfa(nsurfs)='S3'
328  ipkonfa(nsurfs)=ifree
329  do m=1,nopem
330  if(nope.eq.4) then
331  konfa(ifree+m)=konl(ifacet(m,jfacem))
332  else
333  konfa(ifree+m)=konl(ifacew1(m,jfacem))
334  endif
335  enddo
336  ifree=ifree+nopem
337  endif
338  indexe=nodface(5,indexe)
339  if(indexe.eq.0) exit
340  enddo
341  enddo
342 !
343  return
subroutine isortii(ix, iy, n, kflag)
Definition: isortii.f:6
subroutine nodes(inpc, textpart, co, nk, nk_, set, istartset, iendset, ialset, nset, nset_, nalset, nalset_, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: nodes.f:22
Hosted by OpenAircraft.com, (Michigan UAV, LLC)