-
Notifications
You must be signed in to change notification settings - Fork 21
/
manualExplicitIntegrator.py
216 lines (161 loc) · 8.53 KB
/
manualExplicitIntegrator.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
# This is an EXUDYN example
#
# Details: ANCF Cable2D cantilever bent with manual explicit integrator
#
# Author: Johannes Gerstmayr
# Date: 2020-01-08
#
# 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.
#
#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
import exudyn as exu
from exudyn.itemInterface import *
useGraphics = True #without test
#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
#you can erase the following lines and all exudynTestGlobals related operations if this is not intended to be used as TestModel:
try: #only if called from test suite
from modelUnitTests import exudynTestGlobals #for globally storing test results
useGraphics = exudynTestGlobals.useGraphics
except:
class ExudynTestGlobals:
pass
exudynTestGlobals = ExudynTestGlobals()
#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
SC = exu.SystemContainer()
mbs = SC.AddSystem()
exu.Print("\n\n++++++++++++++++++++++++++\nStart EXUDYN version "+exu.GetVersionString()+"\n")
#background
rect = [-2,-2,2,2] #xmin,ymin,xmax,ymax
background0 = {'type':'Line', 'color':[0.1,0.1,0.8,1], 'data':[rect[0],rect[1],0, rect[2],rect[1],0, rect[2],rect[3],0, rect[0],rect[3],0, rect[0],rect[1],0]} #background
background1 = {'type':'Circle', 'radius': 0.1, 'position': [-1.5,0,0]}
background2 = {'type':'Text', 'position': [-1,-1,0], 'text':'Example with text\nin two lines:.=!'} #background
oGround=mbs.AddObject(ObjectGround(referencePosition= [0,0,0], visualization=VObjectGround(graphicsData= [background0, background1, background2])))
#+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
#cable:
mypi = 3.141592653589793
L=2. # length of ANCF element in m
#L=mypi # length of ANCF element in m
E=2.07e11*1e-5 # Young's modulus of ANCF element in N/m^2
rho=7800 # density of ANCF element in kg/m^3
b=0.1 # width of rectangular ANCF element in m
h=0.1 # height of rectangular ANCF element in m
A=b*h # cross sectional area of ANCF element in m^2
I=b*h**3/12 # second moment of area of ANCF element in m^4
f=3*E*I/L**2 # tip load applied to ANCF element in N
exu.Print("load f="+str(f))
exu.Print("EI="+str(E*I))
nGround = mbs.AddNode(NodePointGround(referenceCoordinates=[0,0,0])) #ground node for coordinate constraint
mGround = mbs.AddMarker(MarkerNodeCoordinate(nodeNumber = nGround, coordinate=0)) #Ground node ==> no action
cableList=[]
nc0 = mbs.AddNode(Point2DS1(referenceCoordinates=[0,0,1,0]))
nElements = 4 #for tests use nElements = 4
lElem = L / nElements
nLast = 0
for i in range(nElements):
nLast = mbs.AddNode(Point2DS1(referenceCoordinates=[lElem*(i+1),0,1,0]))
elem=mbs.AddObject(Cable2D(physicsLength=lElem, physicsMassPerLength=rho*A,
physicsBendingStiffness=E*I, physicsAxialStiffness=E*A*0.1,
nodeNumbers=[int(nc0)+i,int(nc0)+i+1], useReducedOrderIntegration=True))
cableList+=[elem]
mANCF0 = mbs.AddMarker(MarkerNodeCoordinate(nodeNumber = nc0, coordinate=0))
mANCF1 = mbs.AddMarker(MarkerNodeCoordinate(nodeNumber = nc0, coordinate=1))
mANCF2 = mbs.AddMarker(MarkerNodeCoordinate(nodeNumber = nc0, coordinate=3))
#mbs.AddObject(CoordinateConstraint(markerNumbers=[mGround,mANCF0]))
#mbs.AddObject(CoordinateConstraint(markerNumbers=[mGround,mANCF1]))
#mbs.AddObject(CoordinateConstraint(markerNumbers=[mGround,mANCF2]))
mANCFLast = mbs.AddMarker(MarkerNodePosition(nodeNumber=nLast)) #force
mbs.AddLoad(Force(markerNumber = mANCFLast, loadVector = [0, -1000, 0])) #will be changed in load steps
#mANCFrigid = mbs.AddMarker(MarkerBodyRigid(bodyNumber=elem, localPosition=[lElem,0,0])) #local position L = beam tip
#mbs.AddLoad(Torque(markerNumber = mANCFrigid, loadVector = [0, 0, E*I*0.25*mypi]))
#mANCFnode = mbs.AddMarker(MarkerNodeRigid(nodeNumber=nLast)) #local position L = beam tip
#mbs.AddLoad(Torque(markerNumber = mANCFnode, loadVector = [0, 0, 0.4*E*I*0.25*mypi]))
#mbs.AddLoad(Force(markerNumber = mANCFnode, loadVector = [0, 0.4*E*I*0.25*mypi,0]))
mbs.Assemble()
#exu.Print(mbs)
simulationSettings = exu.SimulationSettings() #takes currently set values or default values
#SC.visualizationSettings.bodies.showNumbers = False
SC.visualizationSettings.nodes.defaultSize = 0.025
dSize=0.01
SC.visualizationSettings.bodies.defaultSize = [dSize, dSize, dSize]
#simulationSettings.staticSolver.newton.numericalDifferentiation.relativeEpsilon = 1e-9
simulationSettings.staticSolver.verboseMode = 1
simulationSettings.staticSolver.verboseModeFile = 2
simulationSettings.solutionSettings.solverInformationFileName = 'solution/solverInformation.txt'
#simulationSettings.staticSolver.newton.absoluteTolerance = 1e-8
simulationSettings.staticSolver.newton.relativeTolerance = 1e-6 #1e-5 works for 64 elements
simulationSettings.staticSolver.newton.maxIterations = 20 #50 for bending into circle
if useGraphics: #only start graphics once, but after background is set
exu.StartRenderer()
simulationSettings.staticSolver.numberOfLoadSteps = 10
simulationSettings.staticSolver.adaptiveStep = True
import numpy as np
testRefVal = 0
#compute eigenvalues manually:
calcEig = True
if calcEig:
from scipy.linalg import solve, eigh, eig #eigh for symmetric matrices, positive definite
staticSolver = exu.MainSolverStatic()
#staticSolver.SolveSystem(mbs, simulationSettings)
staticSolver.InitializeSolver(mbs, simulationSettings)
staticSolver.ComputeMassMatrix(mbs)
m = staticSolver.GetSystemMassMatrix()
#exu.Print("m =",m)
staticSolver.ComputeJacobianODE2RHS(mbs, scalarFactor_ODE2=-1, scalarFactor_ODE2_t=0)
staticSolver.ComputeJacobianAE(mbs)
K = staticSolver.GetSystemJacobian()
#exu.Print("K =",K)
nODE2 = staticSolver.GetODE2size()
K2 = K[0:nODE2,0:nODE2]
[eigvals, eigvecs] = eigh(K2, m) #this gives omega^2 ... squared eigen frequencies (rad/s)
ev = np.sort(a=abs(eigvals))
#exu.Print("ev =",ev)
if (len(ev) >= 7):
f6 = np.sqrt(abs(ev[6]))/(2*np.pi)
exu.Print("ev=", f6)
testRefVal += f6 #first bending eigenmode
staticSolver.FinalizeSolver(mbs, simulationSettings)
#++++++++++++++++++++++++++++++++++++++++++++++++++
#TEST
def UserFunctionInitializeStep(mainSolver, mainSys, sims):
#exu.Print("t=", mainSolver.it.currentTime)
mainSolver.UpdateCurrentTime(mainSys, sims)
mainSys.systemData.SetTime(mainSolver.it.currentTime);
return True
#test for explicit integrator:
def UserFunctionNewton(mainSolver, mainSys, sims):
nODE2 = mainSolver.GetODE2size()
nAE = mainSolver.GetAEsize()
#nSys = nODE2+nAE
#print("u=", mainSys.systemData.GetODE2Coordinates())
dynamicSolver.ComputeODE2RHS(mbs)
res = dynamicSolver.GetSystemResidual()
Fode2 = res[0:nODE2]
#print("res=", Fode2)
dynamicSolver.ComputeMassMatrix(mbs)
M = dynamicSolver.GetSystemMassMatrix()
a = np.linalg.solve(M,Fode2) #acceleration
h = dynamicSolver.it.currentStepSize
u0 = mainSys.systemData.GetODE2Coordinates()
v0 = mainSys.systemData.GetODE2Coordinates_t()
mainSys.systemData.SetODE2Coordinates(u0+h*v0)
mainSys.systemData.SetODE2Coordinates_t(v0+h*a)
return True
dynamicSolver = exu.MainSolverImplicitSecondOrder()
simulationSettings.timeIntegration.numberOfSteps = 5000 #1000 steps for test suite/error
simulationSettings.timeIntegration.endTime = 0.05 #1s for test suite / error
simulationSettings.timeIntegration.generalizedAlpha.spectralRadius = 0.5
#simulationSettings.displayComputationTime = True
simulationSettings.timeIntegration.verboseMode = 1
#dynamicSolver.SetUserFunctionInitializeStep(mbs, UserFunctionInitializeStep)
dynamicSolver.SetUserFunctionNewton(mbs, UserFunctionNewton)
dynamicSolver.SolveSystem(mbs, simulationSettings)
#mbs.SolveDynamic(simulationSettings)
uy=mbs.GetNodeOutput(nLast,exu.OutputVariableType.Position)[1] #y-coordinate of tip
exu.Print("uy=", uy)
exu.Print("testResult=", testRefVal + uy)
exudynTestGlobals.testError = testRefVal + uy - (2.280183538481952-0.2204849087896498) #2020-01-16: 2.280183538481952-0.2204849087896498
exudynTestGlobals.testResult = testRefVal + uy
if useGraphics: #only start graphics once, but after background is set
SC.WaitForRenderEngineStopFlag()
exu.StopRenderer() #safely close rendering window!