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

Go to the source code of this file.

Functions/Subroutines

subroutine normalsonsurface_se (ipkon, kon, lakon, extnor, co, nk, ipoface, nodface, nactdof, mi, nodedesiinv, noregion)
 

Function/Subroutine Documentation

◆ normalsonsurface_se()

subroutine normalsonsurface_se ( integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
real*8, dimension(3,*)  extnor,
real*8, dimension(3,*)  co,
integer  nk,
integer, dimension(*)  ipoface,
integer, dimension(5,*)  nodface,
integer, dimension(0:mi(2),*)  nactdof,
integer, dimension(*)  mi,
integer, dimension(*)  nodedesiinv,
integer  noregion 
)
21 !
22  implicit none
23 !
24  character*8 lakon(*)
25 !
26  integer j,nelemm,jfacem,indexe,ipkon(*),kon(*),nopem,node,
27  & ifaceq(8,6),ifacet(6,4),ifacew1(4,5),ifacew2(8,5),
28  & konl(26),ipoface(*),nodface(5,*),mi(*),nodedesiinv(*),
29  & nactdof(0:mi(2),*),nopesurf(9),nnodes,noregion,nope,
30  & nopedesi,l,m,iflag,k,nk
31 !
32  real*8 extnor(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 ! calculation of the normal to the external surface;
90 ! each external face contributes its normal
91 ! to each node belonging to the face providing that
92 ! 1) the node is not a design variable OR
93 ! 2) the node is a design variable and belongs to a design face
94 ! A design face is an external face for which more than nopedesi
95 ! nodes are design variables
96 !
97  do j=1,nk
98 !
99  if(ipoface(j).eq.0) cycle
100  indexe=ipoface(j)
101 !
102  do
103 !
104  nelemm=nodface(3,indexe)
105  jfacem=nodface(4,indexe)
106 !
107 ! nopem: # of nodes in the surface
108 ! nope: # of nodes in the element
109 !
110  if(lakon(nelemm)(4:4).eq.'8') then
111  nopem=4
112  nope=8
113  nopedesi=3
114  elseif(lakon(nelemm)(4:5).eq.'20') then
115  nopem=8
116  nope=20
117  nopedesi=5
118  elseif(lakon(nelemm)(4:5).eq.'10') then
119  nopem=6
120  nope=10
121  nopedesi=3
122  elseif(lakon(nelemm)(4:4).eq.'4') then
123  nopem=3
124  nope=4
125  nopedesi=3
126 !
127 ! treatment of wedge faces
128 !
129  elseif(lakon(nelemm)(4:4).eq.'6') then
130  nope=6
131  if(jfacem.le.2) then
132  nopem=3
133  else
134  nopem=4
135  endif
136  nopedesi=3
137  elseif(lakon(nelemm)(4:5).eq.'15') then
138  nope=15
139  if(jfacem.le.2) then
140  nopem=6
141  else
142  nopem=8
143  endif
144  nopedesi=3
145  endif
146  if(noregion.eq.1) nopedesi=0
147 !
148 ! actual position of the nodes belonging to the
149 ! master surface
150 !
151  do k=1,nope
152  konl(k)=kon(ipkon(nelemm)+k)
153  enddo
154 !
155  if((nope.eq.20).or.(nope.eq.8)) then
156  do m=1,nopem
157  nopesurf(m)=konl(ifaceq(m,jfacem))
158  do k=1,3
159  xl2m(k,m)=co(k,nopesurf(m))
160  enddo
161  enddo
162  elseif((nope.eq.10).or.(nope.eq.4))
163  & then
164  do m=1,nopem
165  nopesurf(m)=konl(ifacet(m,jfacem))
166  do k=1,3
167  xl2m(k,m)=co(k,nopesurf(m))
168  enddo
169  enddo
170  elseif(nope.eq.15) then
171  do m=1,nopem
172  nopesurf(m)=konl(ifacew2(m,jfacem))
173  do k=1,3
174  xl2m(k,m)=co(k,nopesurf(m))
175  enddo
176  enddo
177  else
178  do m=1,nopem
179  nopesurf(m)=konl(ifacew1(m,jfacem))
180  do k=1,3
181  xl2m(k,m)=co(k,nopesurf(m))
182  enddo
183  enddo
184  endif
185 !
186 ! sum up how many designvariables are on that surface
187 !
188  nnodes=0
189  do m=1,nopem
190  if(nodedesiinv(nopesurf(m)).eq.1) then
191  nnodes=nnodes+1
192  endif
193  enddo
194 !
195 ! calculate the normal vector in the nodes belonging to the master surface
196 !
197  if(nopem.eq.8) then
198  do m=1,nopem
199  xi=xquad(1,m)
200  et=xquad(2,m)
201  call shape8q(xi,et,xl2m,xsj2,xs2,shp2,iflag)
202  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
203  & + xsj2(3)*xsj2(3))
204  xsj2(1)=xsj2(1)/dd
205  xsj2(2)=xsj2(2)/dd
206  xsj2(3)=xsj2(3)/dd
207 !
208  if(nope.eq.20) then
209  node=konl(ifaceq(m,jfacem))
210  elseif(nope.eq.15) then
211  node=konl(ifacew2(m,jfacem))
212  endif
213  if((nodedesiinv(node).eq.0).or.
214  & ((nodedesiinv(node).eq.1).and.
215  & (nnodes.ge.nopedesi))) then
216  extnor(1,node)=extnor(1,node)
217  & +xsj2(1)
218  extnor(2,node)=extnor(2,node)
219  & +xsj2(2)
220  extnor(3,node)=extnor(3,node)
221  & +xsj2(3)
222  endif
223  enddo
224  elseif(nopem.eq.4) then
225  do m=1,nopem
226  xi=xquad(1,m)
227  et=xquad(2,m)
228  call shape4q(xi,et,xl2m,xsj2,xs2,shp2,iflag)
229  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
230  & + xsj2(3)*xsj2(3))
231  xsj2(1)=xsj2(1)/dd
232  xsj2(2)=xsj2(2)/dd
233  xsj2(3)=xsj2(3)/dd
234 !
235  if(nope.eq.8) then
236  node=konl(ifaceq(m,jfacem))
237  elseif(nope.eq.6) then
238  node=konl(ifacew1(m,jfacem))
239  endif
240  if((nodedesiinv(node).eq.0).or.
241  & ((nodedesiinv(node).eq.1).and.
242  & (nnodes.ge.nopedesi))) then
243  extnor(1,node)=extnor(1,node)
244  & +xsj2(1)
245  extnor(2,node)=extnor(2,node)
246  & +xsj2(2)
247  extnor(3,node)=extnor(3,node)
248  & +xsj2(3)
249  endif
250  enddo
251  elseif(nopem.eq.6) then
252  do m=1,nopem
253  xi=xtri(1,m)
254  et=xtri(2,m)
255  call shape6tri(xi,et,xl2m,xsj2,xs2,shp2,iflag)
256  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
257  & + xsj2(3)*xsj2(3))
258  xsj2(1)=xsj2(1)/dd
259  xsj2(2)=xsj2(2)/dd
260  xsj2(3)=xsj2(3)/dd
261 !
262  if(nope.eq.10) then
263  node=konl(ifacet(m,jfacem))
264  elseif(nope.eq.15) then
265  node=konl(ifacew2(m,jfacem))
266  endif
267  if((nodedesiinv(node).eq.0).or.
268  & ((nodedesiinv(node).eq.1).and.
269  & (nnodes.ge.nopedesi))) then
270  extnor(1,node)=extnor(1,node)
271  & +xsj2(1)
272  extnor(2,node)=extnor(2,node)
273  & +xsj2(2)
274  extnor(3,node)=extnor(3,node)
275  & +xsj2(3)
276  endif
277  enddo
278  else
279  do m=1,nopem
280  xi=xtri(1,m)
281  et=xtri(2,m)
282  call shape3tri(xi,et,xl2m,xsj2,xs2,shp2,iflag)
283  dd=dsqrt(xsj2(1)*xsj2(1) + xsj2(2)*xsj2(2)
284  & + xsj2(3)*xsj2(3))
285  xsj2(1)=xsj2(1)/dd
286  xsj2(2)=xsj2(2)/dd
287  xsj2(3)=xsj2(3)/dd
288 !
289  if(nope.eq.6) then
290  node=konl(ifacew1(m,jfacem))
291  elseif(nope.eq.4) then
292  node=konl(ifacet(m,jfacem))
293  endif
294  if((nodedesiinv(node).eq.0).or.
295  & ((nodedesiinv(node).eq.1).and.
296  & (nnodes.ge.nopedesi))) then
297  extnor(1,node)=extnor(1,node)
298  & +xsj2(1)
299  extnor(2,node)=extnor(2,node)
300  & +xsj2(2)
301  extnor(3,node)=extnor(3,node)
302  & +xsj2(3)
303  endif
304  enddo
305  endif
306 !
307  indexe=nodface(5,indexe)
308  if(indexe.eq.0) exit
309 !
310  enddo
311  enddo
312 !
313 ! not considering the directions with fixed displacements
314 !
315  do l=1,nk
316  do m=1,3
317  if(nactdof(m,l).le.0) then
318  extnor(m,l)=0.d0
319  endif
320  enddo
321  enddo
322 !
323 ! normalizing the normals
324 !
325  do l=1,nk
326  dd=dsqrt(extnor(1,l)**2+extnor(2,l)**2+
327  & extnor(3,l)**2)
328  if(dd.gt.0.d0) then
329  do m=1,3
330  extnor(m,l)=extnor(m,l)/dd
331  enddo
332  endif
333  enddo
334 !
335  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 shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)