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