ANCFBeamTest.py
You can view and download this file on Github: ANCFBeamTest.py
1#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2# This is an EXUDYN example
3#
4# Details: Test models for ANCFBeam (2-node 3D shear deformable ANCF beam element
5# with 2 cross section slope vectors);
6# test models: cantilever beam with tip force and torque
7#
8# Author: Johannes Gerstmayr
9# Date: 2023-04-05
10#
11# Copyright:This file is part of Exudyn. Exudyn is free software. You can redistribute it and/or modify it under the terms of the Exudyn license. See 'LICENSE.txt' for more details.
12#
13#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
14
15import exudyn as exu
16from exudyn.utilities import * #includes itemInterface and rigidBodyUtilities
17import exudyn.graphics as graphics #only import if it does not conflict
18
19import numpy as np
20
21useGraphics = True #without test
22#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23#you can erase the following lines and all exudynTestGlobals related operations if this is not intended to be used as TestModel:
24try: #only if called from test suite
25 from modelUnitTests import exudynTestGlobals #for globally storing test results
26 useGraphics = exudynTestGlobals.useGraphics
27except:
28 class ExudynTestGlobals:
29 pass
30 exudynTestGlobals = ExudynTestGlobals()
31#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32
33SC = exu.SystemContainer()
34mbs = SC.AddSystem()
35
36compute2D = False
37compute3D = True
38
39#test examples
40#2011 MUBO, Nachbagauer Pechstein Irschik Gerstmayr (2D)
41#2013 CND, Nachbagauer Gruber Gerstmayr (static, 3D); "Structural and Continuum Mechanics Approaches for a 3D Shear Deformable ANCF Beam Finite Element: Application to Static and Linearized Dynamic Examples"
42### not yet: 2013 CND, Nachbagauer Gerstmayr (dynamic, 3D)
43cases = ['CantileverLinear2011', 'Cantilever2011', 'GeneralBending2013', 'PrincetonBeamF2', 'PrincetonBeamF3', 'Eigenmodes2013']
44nElementsList = [1,2,4,8,16,32,64,128,256,512,1024]
45# nElementsList = [8,32, 128]
46nElements = 8
47
48betaList = [0,15,30,45,60,75,90]
49betaDegree = 45
50
51caseList = [0,1,2,3,4]
52case=4
53
54useGraphics = False
55verbose = False
56
57bodyFixedLoad = False
58testErrorSum = 0
59printCase = True
60
61#for nElements in nElementsList:
62#for betaDegree in betaList:
63for case in caseList:
64 if printCase:
65 printCase=False
66 exu.Print('case=', case, cases[case])
67 mbs.Reset()
68
69 computeEigenmodes = False
70 csFact = 1
71 sectionData = exu.BeamSection()
72 fTip = 0
73 MxTip = 0
74 MyTip = 0
75
76 ks1=1 #shear correction, torsion
77 ks2=1 #shear correction, bending
78 ks3=1 #shear correction, bending
79 ff=1 #drawing factor
80
81 #define beam parameters and loads for different cases
82 if case == 0 or case == 1:
83 caseName = cases[case]
84
85 L = 2 #length of beam
86 w = 0.1 #width of beam
87 h = 0.5 #height Y
88
89 fTip = 5e5*h**3
90 if case == 1:
91 fTip *= 1000
92
93 Em = 2.07e11
94 rho = 1e2
95
96 A=h*w
97 nu = 0.3 # Poisson ratio
98 ks2= 10*(1+nu)/(12+11*nu)
99 ks3=ks2
100
101 elif case == 2:
102 L = 2 #length of beam
103 h = 0.2 #height Y
104 w = 0.4 #width Z of beam
105 Em = 2.07e11
106 rho = 1e2
107
108 A=h*w
109
110 nu = 0.3 # Poisson ratio
111 ks1= 0.5768 #torsion correction factor if J=Jyy+Jzz
112 ks2= 0.8331
113 ks3= 0.7961
114
115 MxTip = 0.5e6
116 MyTip = 2e6
117
118 csFact = 10
119 elif case == 3 or case == 4: #Princeton beam example
120 L = 0.508 #length of beam
121 h = 12.3777e-3 #height Y; 12.3777e-3 with Obrezkov's paper
122 w = 3.2024e-3 #width Z of beam
123 Em = 71.7e9
124 ks1=0.198
125 nu = 0.31
126
127 ks2=1
128 ks3=1
129 # ks2=0.9
130 # ks3=0.9
131
132
133 rho = 1e2 #unused
134 A=h*w
135
136 MxTip = 0
137 MyTip = 0
138 if case == 3:
139 fTip = 8.896 #F2
140 elif case == 4:
141 fTip = 13.345 #F3
142 #if kk==0: exu.Print('load=', fTip)
143
144 beta = betaDegree/180*pi #beta=0 => negative y-axis
145 bodyFixedLoad = False
146
147 csFact = 10
148
149 Gm = Em/(2*(1+nu)) # Shear modulus
150
151 #compute sectionData
152
153 # Cross-section properties
154 Iyy = h*w**3/12 # Second moment of area of the beam cross-section
155 Izz = w*h**3/12 # Second moment of area of the beam cross-section
156 J = (Iyy+Izz) # approximation; Polar moment of area of the beam cross-section
157
158 sectionData.stiffnessMatrix = np.diag([Em*A, Gm*A*ks2, Gm*A*ks3, Gm*J*ks1, Em*Iyy, Em*Izz])
159
160
161 rhoA = rho*A
162
163 if False:
164 #linear solution:
165 uzTip = fTip*L**3/(3*Em*Iyy)
166 exu.Print('uz linear=',uzTip)
167 uyTip = fTip*L**3/(3*Em*Izz)
168 exu.Print('uy linear=',uyTip)
169
170 sectionData.inertia= rho*J*np.eye(3)
171 sectionData.massPerLength = rhoA
172
173 sectionGeometry = exu.BeamSectionGeometry()
174
175 #points, in positive rotation sense viewing in x-direction, points in [Y,Z]-plane
176 #points do not need to be closed!
177 lp = exu.Vector2DList()
178 if True:
179 lp.Append([h*ff,-w*ff])
180 lp.Append([h*ff,w*ff])
181 lp.Append([-h*ff,w*ff])
182 lp.Append([-h*ff,-w*ff])
183
184 sectionGeometry.polygonalPoints = lp
185 #exu.Print('HERE\n',sectionGeometry.polygonalPoints)
186 nGround = mbs.AddNode(NodePointGround(referenceCoordinates=[0,0,0])) #ground node for coordinate constraint
187 mnGround = mbs.AddMarker(MarkerNodeCoordinate(nodeNumber=nGround, coordinate=0))
188
189
190 eY=[0,1,0]
191 eZ=[0,0,1]
192 lElem = L/nElements
193 if compute3D:
194 initialRotations = eY+eZ
195 #create beam nodes and elements
196 n0 = mbs.AddNode(NodePointSlope23(referenceCoordinates=[0,0,0]+initialRotations))
197 nInit = n0
198 for k in range(nElements):
199 n1 = mbs.AddNode(NodePointSlope23(referenceCoordinates=[lElem*(k+1),0,0]+initialRotations))
200
201 oBeam = mbs.AddObject(ObjectANCFBeam(nodeNumbers=[n0,n1], physicsLength = lElem,
202 sectionData = sectionData,
203 crossSectionPenaltyFactor = [csFact,csFact,csFact],
204 visualization=VANCFBeam(sectionGeometry=sectionGeometry)))
205 n0 = n1
206
207
208 mTip = mbs.AddMarker(MarkerNodeRigid(nodeNumber = n1))
209 if fTip != 0:
210 if case < 3:
211 mbs.AddLoad(Force(markerNumber=mTip, loadVector = [0,fTip,0], bodyFixed = bodyFixedLoad))
212 elif case >= 3:
213 mbs.AddLoad(Force(markerNumber=mTip, loadVector = [0,-fTip*cos(beta),fTip*sin(beta)], bodyFixed = bodyFixedLoad))
214
215 if MxTip != 0 or MyTip != 0:
216 mbs.AddLoad(Torque(markerNumber=mTip, loadVector = [MxTip, MyTip,0]))#, bodyFixed = True ))
217
218 for i in range(9):
219 #if i != 4 and i != 8: #exclude constraining the slope lengths
220 if True:
221 nm0 = mbs.AddMarker(MarkerNodeCoordinate(nodeNumber=nInit, coordinate=i))
222 mbs.AddObject(CoordinateConstraint(markerNumbers=[mnGround, nm0]))
223
224
225 # exu.Print(mbs)
226 mbs.Assemble()
227
228 tEnd = 100 #end time of simulation
229 stepSize = 0.5*0.01*0.1 #step size; leads to 1000 steps
230
231 simulationSettings = exu.SimulationSettings()
232 simulationSettings.solutionSettings.solutionWritePeriod = 2e-2 #output interval general
233 simulationSettings.solutionSettings.sensorsWritePeriod = 1e-1 #output interval of sensors
234 simulationSettings.timeIntegration.numberOfSteps = int(tEnd/stepSize) #must be integer
235 simulationSettings.timeIntegration.endTime = tEnd
236 #simulationSettings.solutionSettings.solutionInformation = "This is the info\nNew line\n and another new line \n"
237 simulationSettings.timeIntegration.generalizedAlpha.spectralRadius = 0.5
238 #simulationSettings.timeIntegration.simulateInRealtime=True
239 #simulationSettings.timeIntegration.realtimeFactor=0.1
240
241 simulationSettings.timeIntegration.verboseMode = verbose
242 simulationSettings.staticSolver.verboseMode = verbose
243
244 simulationSettings.timeIntegration.newton.useModifiedNewton = True
245
246 simulationSettings.timeIntegration.newton.numericalDifferentiation.relativeEpsilon = 1e-4
247 simulationSettings.timeIntegration.newton.relativeTolerance = 1e-6
248
249 simulationSettings.linearSolverType = exu.LinearSolverType.EigenSparse
250
251 if nElements > 32 and case==0: #change tolerance, because otherwise no convergence
252 simulationSettings.staticSolver.newton.relativeTolerance = 1e-6
253 if case == 1: #tolerance changed from 1e-8 to 5e-10 to achieve values of paper (1024 has difference at last digit in paper)
254 simulationSettings.staticSolver.newton.relativeTolerance = 0.5e-9
255
256
257 simulationSettings.staticSolver.numberOfLoadSteps = 5
258
259
260 #add some drawing parameters for this example
261 SC.visualizationSettings.nodes.drawNodesAsPoint=False
262 SC.visualizationSettings.nodes.defaultSize=0.01
263
264 SC.visualizationSettings.bodies.beams.axialTiling = 50
265 SC.visualizationSettings.general.drawWorldBasis = True
266 SC.visualizationSettings.general.worldBasisSize = 0.1
267 SC.visualizationSettings.openGL.multiSampling = 4
268
269
270 # [M, K, D] = exu.solver.ComputeLinearizedSystem(mbs, simulationSettings, useSparseSolver=True)
271 # exu.Print('M=',M.round(1))
272
273 if useGraphics:
274 exu.StartRenderer()
275 mbs.WaitForUserToContinue()
276
277 # if computeEigenmodes:
278 # nModes = 3*(1+int(compute3D))
279 # nRigidModes = 3*(1+int(compute3D))
280 # if compute2D:
281 # constrainedCoordinates=[0,1,mbs.systemData.ODE2Size()-2]
282 # else:
283 # constrainedCoordinates=[0,1,2,5,mbs.systemData.ODE2Size()-8,mbs.systemData.ODE2Size()-7]
284
285 # # constrainedCoordinates=[]
286
287 # compeig=mbs.ComputeODE2Eigenvalues(simulationSettings, useSparseSolver=False,
288 # numberOfEigenvalues= nRigidModes+nModes,
289 # constrainedCoordinates=constrainedCoordinates,
290 # convert2Frequencies= False)
291
292 # exu.Print('eigvalues=',np.sqrt(compeig[0][nRigidModes:]))
293
294 # if False: #show modes:
295 # for i in range(nModes):
296 # iMode = nRigidModes+i
297 # mbs.systemData.SetODE2Coordinates(5*compeig[1][:,iMode], exudyn.ConfigurationType.Visualization)
298 # mbs.systemData.SetTime(np.sqrt(compeig[0][iMode]), exudyn.ConfigurationType.Visualization)
299 # mbs.SendRedrawSignal()
300
301 # mbs.WaitForUserToContinue()
302
303 # else:
304 mbs.SolveStatic(simulationSettings)
305 # mbs.SolveDynamic(simulationSettings)
306 #mbs.SolveDynamic(simulationSettings, solverType = exu.DynamicSolverType.RK44)
307
308
309 if useGraphics:
310 SC.WaitForRenderEngineStopFlag()
311 exu.StopRenderer() #safely close rendering window!
312
313 ##evaluate final (=current) output values
314 uTip = mbs.GetNodeOutput(n1, exu.OutputVariableType.Displacement)
315
316 errorFact = 1
317 if case != 1:
318 errorFact *= 100
319
320 testErrorSum += np.linalg.norm(uTip)
321
322
323
324 if case < 2:
325 pTip = mbs.GetNodeOutput(n1, exu.OutputVariableType.Position)
326 exu.Print('ne=',nElements, ', ux=',L-pTip[0], ', uy=',pTip[1])
327 elif case == 2:
328 rotTip = mbs.GetNodeOutput(n1, exu.OutputVariableType.Rotation)
329 exu.Print('ne=',nElements, ', u=',list(uTip))
330 # exu.Print('ne=',nElements, ', rot=',rotTip)
331 elif case == 3 or case == 4:
332 exu.Print('ne=', nElements, ', beta=', round(beta*180/pi,1), ', u=',uTip.round(7))
333
334
335exu.Print('Solution of ANCFBeam3Dtest=', testErrorSum)
336exudynTestGlobals.testError = testErrorSum - (1.010486312300459 )
337exudynTestGlobals.testResult = testErrorSum
338
339
340
341# case= 0/CantileverLinear2011
342#NachbagauerPechsteinIrschikGerstmayrMUBO2011 (2D):
343# ne=1, 9.12273046e–8, 6.16666566e–4, 0.000193
344# ne=2, 1.61293091e–7, 7.61594059e–4, 4.831e–5
345# ne=4, 1.81763233e–7, 7.97825954e–4, 1.208e–5
346# ne=256, 1.88847418e–7, 8.09900305e–4, 2.945e–9
347#Exudyn: ksFact=1
348# ne= 1 , ux= 9.122730637578513e-08 , uy= 0.0006166665660910789
349# ne= 2 , ux= 1.612930911054633e-07 , uy= 0.0007615940599560586
350# ne= 4 , ux= 1.8176323512975046e-07 , uy= 0.0007978259537503566
351# ne= 8 , ux= 1.8706537496804287e-07 , uy= 0.0008068839288072378
352# ne= 16 , ux= 1.8840244964124508e-07 , uy= 0.0008091484226773518
353# ne= 32 , ux= 1.887374359021976e-07 , uy= 0.0008097145461515286
354# ne= 64 , ux= 1.888212299849812e-07 , uy= 0.0008098560770202866
355# ne= 128 , ux= 1.8884218011550047e-07 , uy= 0.000809891459736643
356# ne= 256 , ux= 1.8884741770364144e-07 , uy= 0.0008099003054122335
357
358
359# case= 1/Cantilever2011
360#NachbagauerPechsteinIrschikGerstmayrMUBO2011 (2D):
361# ne=1, 0.07140274, 0.54225823, 0.168310
362# ne=2, 0.12379212, 0.65687111, 0.053697
363# ne=4, 0.14346767, 0.69593561, 0.014633
364# ne=1024, 0.15097103, 0.71056837, 2.280e–7
365
366#Exudyn: ksFact=1
367# ne= 1 , ux= 0.07140273975041422 , uy= 0.5422582285449739
368# ne= 2 , ux= 0.12379212054619537 , uy= 0.6568711099777776
369# ne= 4 , ux= 0.14346766617229956 , uy= 0.695935613449867
370# ne= 8 , ux= 0.14904162148449163 , uy= 0.7068152604035266
371# ne= 16 , ux= 0.15048521526298897 , uy= 0.709623891842095
372# ne= 32 , ux= 0.15084943688011565 , uy= 0.7103320154655514
373# ne= 64 , ux= 0.15094070328691145 , uy= 0.7105094267817303
374# ne= 128 , ux= 0.15096353326024237 , uy= 0.7105538037895819
375# ne= 256 , ux= 0.15096924149743085 , uy= 0.7105648993600513
376# ne= 512 , ux= 0.15097066651939461 , uy= 0.7105676689547459
377# ne= 1024 , ux= 0.15097102364723924 , uy= 0.7105683631862169
378
379# case = 2:
380#2013 CND, Nachbagauer Gruber Gerstmayr (static, 3D); "Structural and Continuum Mechanics Approaches for a 3D Shear Deformable ANCF Beam Finite Element: Application to Static and Linearized Dynamic Examples"
381#Table 4:
382# SMF
383# 8, 1.0943e-4, 1.8638e-4, 1.8117e-2
384# 32, 1.0943e-4, 1.8625e-4, 1.8117e-2
385# ANSYS
386# 40, 1.0939e-4, 1.8646e-4, 1.8117e-2
387#Exudyn, ksFact=10:
388# ne= 8 , u= [-0.00010900977088157404, -0.0001902100873246334, -0.01811732779800177]
389# ne= 32 , u= [-0.00010941122286522997, -0.00018667435478355072, -0.01811739809277171]
390# ne= 128 , u= [-0.00010943631319815239, -0.000186451835025629, -0.018117402461210096]
391#==> in 2013 paper, element performed slightly better, especially in ux and uy terms
392
393# case = 3:
394#Princeton beam with ANSYS (Leonid Obrezkov / Aki Mikkola / Marko Matikainen et al.,
395# Performance review of locking alleviation methods for continuum ANCF beam elements,
396# Nonlinear Dynamics, Vol. 109, pp. 31–546, May 2022
397# beta=[0 15 30 45 60 75 90];
398if (case==3 or case == 4) and False:
399 # F2=8.896
400 # % ANSYS beam (10-199 el)
401 ANSYSF2y=np.array([1.071417630E-002, 1.061328706E-002, 1.011169630E-002, 8.837226265E-003, 6.604665004E-003, 3.538889001E-003, 0])
402 ANSYSF2z=np.array([0, 4.208232124E-002, 7.939482948E-002, 0.108987937, 0.129887616, 0.142194370, 0.146245978])
403 exu.Print('refsol ANSYS F2=8.896:\n',ANSYSF2y.round(6), '\n', ANSYSF2z.round(6),sep='')
404 # % ANSYS solid (el) (4x12x500) - finer mesh doesn't have much influence see in Size effect file
405 # ANSYS_solid_y=[1.069752828E-002 1.057180106E-002 9.938278402E-003 8.686786771E-003 6.500006282E-003 3.481999513E-003 0];
406 # ANSYS_solid_z=[0 4.101165651E-002 7.696749069E-002 0.105976311 0.127251299 0.139594740 0.143848652];
407
408 # F3=13.345
409 # % ANSYS beam (10-199 el)
410 ANSYSF3y=np.array([1.606423724E-002, 1.645825752E-002, 1.665873206E-002, 1.518618440E-002, 1.157837500E-002, 6.248967384E-003, 0])
411 ANSYSF3z=np.array([0, 6.435812858E-002, 0.117735994, 0.156467239, 0.181861627, 0.196097131, 0.200677707])
412 # % ANSYS solid (el) (4x12x500) - finer mesh see in Size effect file
413 #ANSYS_solid_y=[1.603700622E-002 1.637026068E-002 1.640440775E-002 1.485055210E-002 1.127173264E-002 6.062461977E-003 0])
414 #ANSYS_solid_z=[0 6.270699533E-002 0.113752002 0.153554457 0.179978534 0.192972233 0.197669499])
415 exu.Print('refsol ANSYS F3=13.345:\n',ANSYSF3y.round(6), '\n', ANSYSF3z.round(6),sep='')
416#Exudyn results for Princeton beam:
417#not exactly the same, but around the previous values with HOTINT
418#using 16 elements, csFact=10 (no influence)
419# F2=8.896
420# case= 3, PrincetonBeam
421# ne= 16 , beta= 0.0 , u= [-0.0001352 -0.0107023 0. ]
422# ne= 16 , beta= 15.0 , u= [-0.0022414 -0.0106295 0.0421374]
423# ne= 16 , beta= 30.0 , u= [-0.0076567 -0.0101861 0.0794434]
424# ne= 16 , beta= 45.0 , u= [-0.0143664 -0.0089529 0.1089703]
425# ne= 16 , beta= 60.0 , u= [-0.0204225 -0.0067182 0.1297877]
426# ne= 16 , beta= 75.0 , u= [-0.0245093 -0.0036079 0.1420319]
427# ne= 16 , beta= 90.0 , u= [-0.0259403 -0. 0.1460608]
428
429# F3=13.345
430# case= 4, PrincetonBeam
431# ne= 16 , beta= 0.0 , u= [-0.0003039 -0.0160454 0. ]
432# ne= 16 , beta= 15.0 , u= [-0.005319 -0.0165469 0.064622 ]
433# ne= 16 , beta= 30.0 , u= [-0.0171901 -0.0169316 0.1179818]
434# ne= 16 , beta= 45.0 , u= [-0.0303357 -0.0155488 0.1565214]
435# ne= 16 , beta= 60.0 , u= [-0.0411035 -0.0118996 0.1817173]
436# ne= 16 , beta= 75.0 , u= [-0.0479101 -0.0064334 0.1958343]
437# ne= 16 , beta= 90.0 , u= [-0.0502184 -0. 0.2003738]