-
Notifications
You must be signed in to change notification settings - Fork 299
/
EulerImplicitSolver.cpp
350 lines (292 loc) · 13.6 KB
/
EulerImplicitSolver.cpp
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
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
/******************************************************************************
* SOFA, Simulation Open-Framework Architecture *
* (c) 2006 INRIA, USTL, UJF, CNRS, MGH *
* *
* This program is free software; you can redistribute it and/or modify it *
* under the terms of the GNU Lesser General Public License as published by *
* the Free Software Foundation; either version 2.1 of the License, or (at *
* your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, but WITHOUT *
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or *
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License *
* for more details. *
* *
* You should have received a copy of the GNU Lesser General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*******************************************************************************
* Authors: The SOFA Team and external contributors (see Authors.txt) *
* *
* Contact information: contact@sofa-framework.org *
******************************************************************************/
#include <sofa/component/odesolver/backward/EulerImplicitSolver.h>
#include <sofa/core/visual/VisualParams.h>
#include <sofa/simulation/MechanicalOperations.h>
#include <sofa/simulation/VectorOperations.h>
#include <sofa/helper/AdvancedTimer.h>
#include <sofa/core/ObjectFactory.h>
#include <sofa/core/behavior/MultiMatrix.h>
#include <sofa/helper/ScopedAdvancedTimer.h>
namespace sofa::component::odesolver::backward
{
using core::VecId;
using namespace sofa::defaulttype;
using namespace core::behavior;
EulerImplicitSolver::EulerImplicitSolver()
: d_rayleighStiffness(initData(&d_rayleighStiffness, (SReal)0.0, "rayleighStiffness", "Rayleigh damping coefficient related to stiffness, > 0") )
, d_rayleighMass(initData(&d_rayleighMass, (SReal)0.0, "rayleighMass", "Rayleigh damping coefficient related to mass, > 0"))
, d_velocityDamping(initData(&d_velocityDamping, (SReal)0.0, "vdamping", "Velocity decay coefficient (no decay if null)") )
, d_firstOrder (initData(&d_firstOrder, false, "firstOrder", "Use backward Euler scheme for first order ode system."))
, d_trapezoidalScheme( initData(&d_trapezoidalScheme,false,"trapezoidalScheme","Optional: use the trapezoidal scheme instead of the implicit Euler scheme and get second order accuracy in time") )
, d_solveConstraint(initData(&d_solveConstraint, false, "solveConstraint", "Apply ConstraintSolver (requires a ConstraintSolver in the same node as this solver, disabled by by default for now)") )
, d_threadSafeVisitor(initData(&d_threadSafeVisitor, false, "threadSafeVisitor", "If true, do not use realloc and free visitors in fwdInteractionForceField."))
{
f_rayleighStiffness.setParent(&d_rayleighStiffness);
f_rayleighMass.setParent(&d_rayleighMass);
f_velocityDamping.setParent(&d_velocityDamping);
f_firstOrder.setParent(&d_firstOrder);
f_solveConstraint.setParent(&d_solveConstraint);
}
void EulerImplicitSolver::init()
{
if (!this->getTags().empty())
{
msg_info() << "EulerImplicitSolver: responsible for the following objects with tags " << this->getTags() << " :";
type::vector<core::objectmodel::BaseObject*> objs;
this->getContext()->get<core::objectmodel::BaseObject>(&objs,this->getTags(),sofa::core::objectmodel::BaseContext::SearchDown);
for (const auto* obj : objs)
msg_info() << " " << obj->getClassName() << ' ' << obj->getName();
}
sofa::core::behavior::OdeSolver::init();
}
void EulerImplicitSolver::cleanup()
{
// free the locally created vector x (including eventual external mechanical states linked by an InteractionForceField)
sofa::simulation::common::VectorOperations vop( core::execparams::defaultInstance(), this->getContext() );
vop.v_free(x.id(), !d_threadSafeVisitor.getValue(), true);
}
void EulerImplicitSolver::solve(const core::ExecParams* params, SReal dt, sofa::core::MultiVecCoordId xResult, sofa::core::MultiVecDerivId vResult)
{
#ifdef SOFA_DUMP_VISITOR_INFO
sofa::simulation::Visitor::printNode("SolverVectorAllocation");
#endif
sofa::simulation::common::VectorOperations vop( params, this->getContext() );
sofa::simulation::common::MechanicalOperations mop( params, this->getContext() );
MultiVecCoord pos(&vop, core::VecCoordId::position() );
MultiVecDeriv vel(&vop, core::VecDerivId::velocity() );
MultiVecDeriv f(&vop, core::VecDerivId::force() );
MultiVecDeriv b(&vop, true, core::VecIdProperties{"RHS", GetClass()->className});
MultiVecCoord newPos(&vop, xResult );
MultiVecDeriv newVel(&vop, vResult );
/// inform the constraint parameters about the position and velocity id
mop.cparams.setX(xResult);
mop.cparams.setV(vResult);
// dx is no longer allocated by default (but it will be deleted automatically by the mechanical objects)
MultiVecDeriv dx(&vop, core::VecDerivId::dx());
dx.realloc(&vop, !d_threadSafeVisitor.getValue(), true);
x.realloc(&vop, !d_threadSafeVisitor.getValue(), true, core::VecIdProperties{"solution", GetClass()->className});
#ifdef SOFA_DUMP_VISITOR_INFO
sofa::simulation::Visitor::printCloseNode("SolverVectorAllocation");
#endif
const SReal& h = dt;
const bool firstOrder = d_firstOrder.getValue();
// the only difference for the trapezoidal rule is the factor tr = 0.5 for some usages of h
const bool optTrapezoidal = d_trapezoidalScheme.getValue();
SReal tr;
if (optTrapezoidal)
tr = 0.5;
else
tr = 1.0;
msg_info() << "trapezoidal factor = " << tr;
{
SCOPED_TIMER("ComputeForce");
mop->setImplicit(true); // this solver is implicit
// compute the net forces at the beginning of the time step
mop.computeForce(f);
msg_info() << "EulerImplicitSolver, initial f = " << f;
}
sofa::helper::AdvancedTimer::stepBegin("ComputeRHTerm");
if( firstOrder )
{
b.eq(f);
}
else
{
// new more powerful visitors
// force in the current configuration
b.eq(f,1.0/tr); // b = f0
msg_info() << "EulerImplicitSolver, f = " << f;
// add the change of force due to stiffness + Rayleigh damping
mop.addMBKv(b, -d_rayleighMass.getValue(), 1, h + d_rayleighStiffness.getValue()); // b = f0 + ( rm M + B + (h+rs) K ) v
// integration over a time step
b.teq(h*tr); // b = h(f0 + ( rm M + B + (h+rs) K ) v )
}
msg_info() << "EulerImplicitSolver, b = " << b;
mop.projectResponse(b); // b is projected to the constrained space
msg_info() << "EulerImplicitSolver, projected b = " << b;
sofa::helper::AdvancedTimer::stepNext ("ComputeRHTerm", "MBKBuild");
core::behavior::MultiMatrix<simulation::common::MechanicalOperations> matrix(&mop);
if (firstOrder)
matrix.setSystemMBKMatrix(MechanicalMatrix(1,0,-h*tr)); //MechanicalMatrix::K * (-h*tr) + MechanicalMatrix::M;
else
matrix.setSystemMBKMatrix(MechanicalMatrix(1+ tr * h * d_rayleighMass.getValue(), -tr * h, -tr * h * (h + d_rayleighStiffness.getValue()))); // MechanicalMatrix::K * (-tr*h*(h+d_rayleighStiffness.getValue())) + MechanicalMatrix::B * (-tr*h) + MechanicalMatrix::M * (1+tr*h*d_rayleighMass.getValue());
msg_info() << "EulerImplicitSolver, matrix = " << (MechanicalMatrix::K * (-h * (h + d_rayleighStiffness.getValue())) + MechanicalMatrix::M * (1 + h * d_rayleighMass.getValue())) << " = " << matrix;
msg_info() << "EulerImplicitSolver, Matrix K = " << MechanicalMatrix::K;
#ifdef SOFA_DUMP_VISITOR_INFO
simulation::Visitor::printNode("SystemSolution");
#endif
sofa::helper::AdvancedTimer::stepEnd ("MBKBuild");
{
SCOPED_TIMER("MBKSolve");
matrix.solve(x, b); //Call to ODE resolution: x is the solution of the system}
}
#ifdef SOFA_DUMP_VISITOR_INFO
simulation::Visitor::printCloseNode("SystemSolution");
#endif
// mop.projectResponse(x);
// x is the solution of the system
// apply the solution
const bool solveConstraint = d_solveConstraint.getValue();
#ifndef SOFA_NO_VMULTIOP // unoptimized version
if (solveConstraint)
#endif
{
if (firstOrder)
{
const char* prevStep = "UpdateV";
sofa::helper::AdvancedTimer::stepBegin(prevStep);
#define SOFATIMER_NEXTSTEP(s) { sofa::helper::AdvancedTimer::stepNext(prevStep,s); prevStep=s; }
newVel.eq(x); // vel = x
if (solveConstraint)
{
SOFATIMER_NEXTSTEP("CorrectV");
mop.solveConstraint(newVel,core::ConstraintOrder::VEL);
}
SOFATIMER_NEXTSTEP("UpdateX");
newPos.eq(pos, newVel, h); // pos = pos + h vel
if (solveConstraint)
{
SOFATIMER_NEXTSTEP("CorrectX");
mop.solveConstraint(newPos,core::ConstraintOrder::POS);
}
#undef SOFATIMER_NEXTSTEP
sofa::helper::AdvancedTimer::stepEnd (prevStep);
}
else
{
const char* prevStep = "UpdateV";
sofa::helper::AdvancedTimer::stepBegin(prevStep);
#define SOFATIMER_NEXTSTEP(s) { sofa::helper::AdvancedTimer::stepNext(prevStep,s); prevStep=s; }
// vel = vel + x
newVel.eq(vel, x);
if (solveConstraint)
{
SOFATIMER_NEXTSTEP("CorrectV");
mop.solveConstraint(newVel,core::ConstraintOrder::VEL);
}
SOFATIMER_NEXTSTEP("UpdateX");
// pos = pos + h vel
newPos.eq(pos, newVel, h);
if (solveConstraint)
{
SOFATIMER_NEXTSTEP("CorrectX");
mop.solveConstraint(newPos,core::ConstraintOrder::POS);
}
#undef SOFATIMER_NEXTSTEP
sofa::helper::AdvancedTimer::stepEnd (prevStep);
}
}
#ifndef SOFA_NO_VMULTIOP
else
{
typedef core::behavior::BaseMechanicalState::VMultiOp VMultiOp;
VMultiOp ops;
if (firstOrder)
{
ops.resize(2);
ops[0].first = newVel;
ops[0].second.push_back(std::make_pair(x.id(),1.0));
ops[1].first = newPos;
ops[1].second.push_back(std::make_pair(pos.id(),1.0));
ops[1].second.push_back(std::make_pair(newVel.id(),h));
}
else
{
ops.resize(2);
ops[0].first = newVel;
ops[0].second.push_back(std::make_pair(vel.id(),1.0));
ops[0].second.push_back(std::make_pair(x.id(),1.0));
ops[1].first = newPos;
ops[1].second.push_back(std::make_pair(pos.id(),1.0));
ops[1].second.push_back(std::make_pair(newVel.id(),h));
}
SCOPED_TIMER_VARNAME(updateVAndXTimer, "UpdateVAndX");
vop.v_multiop(ops);
if (solveConstraint)
{
{
SCOPED_TIMER_VARNAME(correctVTimer, "CorrectV");
mop.solveConstraint(newVel,core::ConstraintOrder::VEL);
}
{
SCOPED_TIMER_VARNAME(correctXTimer, "CorrectX");
mop.solveConstraint(newPos,core::ConstraintOrder::POS);
}
}
}
#endif
mop.addSeparateGravity(dt, newVel); // v += dt*g . Used if mass wants to add G separately from the other forces to v
if (d_velocityDamping.getValue() != 0.0)
newVel *= exp(-h * d_velocityDamping.getValue());
if( f_printLog.getValue() )
{
mop.projectPosition(newPos);
mop.projectVelocity(newVel);
mop.propagateX(newPos);
mop.propagateV(newVel);
msg_info() << "EulerImplicitSolver, final x = " << newPos;
msg_info() << "EulerImplicitSolver, final v = " << newVel;
mop.computeForce(f);
msg_info() << "EulerImplicitSolver, final f = " << f;
}
}
SReal EulerImplicitSolver::getPositionIntegrationFactor() const
{
return getPositionIntegrationFactor(getContext()->getDt());
}
SReal EulerImplicitSolver::getIntegrationFactor(int inputDerivative, int outputDerivative) const
{
return getIntegrationFactor(inputDerivative, outputDerivative, getContext()->getDt());
}
SReal EulerImplicitSolver::getIntegrationFactor(int inputDerivative, int outputDerivative, SReal dt) const
{
const SReal matrix[3][3] =
{
{ 1, dt, 0},
{ 0, 1, 0},
{ 0, 0, 0}
};
if (inputDerivative >= 3 || outputDerivative >= 3)
return 0;
else
return matrix[outputDerivative][inputDerivative];
}
SReal EulerImplicitSolver::getSolutionIntegrationFactor(int outputDerivative) const
{
return getSolutionIntegrationFactor(outputDerivative, getContext()->getDt());
}
SReal EulerImplicitSolver::getSolutionIntegrationFactor(int outputDerivative, SReal dt) const
{
const SReal vect[3] = { dt, 1, 1/dt};
if (outputDerivative >= 3)
return 0;
else
return vect[outputDerivative];
}
int EulerImplicitSolverClass = core::RegisterObject("Time integrator using implicit backward Euler scheme")
.add< EulerImplicitSolver >()
.addAlias("EulerImplicit")
.addAlias("ImplicitEulerSolver")
.addAlias("ImplicitEuler")
;
} // namespace sofa::component::odesolver::backward