-
Notifications
You must be signed in to change notification settings - Fork 297
/
SparseLDLSolver.inl
292 lines (243 loc) · 10.9 KB
/
SparseLDLSolver.inl
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
/******************************************************************************
* 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 *
******************************************************************************/
#pragma once
#include <sofa/component/linearsolver/direct/SparseLDLSolver.h>
#include <sofa/core/visual/VisualParams.h>
#include <sofa/core/ObjectFactory.h>
#include <sofa/helper/system/thread/CTime.h>
#include <sofa/core/objectmodel/BaseContext.h>
#include <sofa/core/behavior/LinearSolver.h>
#include <cmath>
#include <fstream>
#include <iomanip> // std::setprecision
#include <string>
#include <sofa/simulation/MainTaskSchedulerFactory.h>
#include <sofa/simulation/ParallelForEach.h>
namespace sofa::component::linearsolver::direct
{
template<class TMatrix, class TVector, class TThreadManager>
SparseLDLSolver<TMatrix,TVector,TThreadManager>::SparseLDLSolver()
: numStep(0)
{}
template <class TMatrix, class TVector, class TThreadManager>
void SparseLDLSolver<TMatrix, TVector, TThreadManager>::init()
{
Inherit::init();
if (this->d_componentState.getValue() != core::objectmodel::ComponentState::Invalid)
{
this->d_componentState.setValue(core::objectmodel::ComponentState::Valid);
}
}
template <class TMatrix, class TVector, class TThreadManager>
void SparseLDLSolver<TMatrix, TVector, TThreadManager>::parse(sofa::core::objectmodel::BaseObjectDescription* arg)
{
Inherit1::parse(arg);
if (!arg->getAttribute("template"))
{
std::string header = this->getClassName();
if (const std::string& name = this->getName(); !name.empty())
{
header.append("(" + name + ")");
}
static const char* blocksType =
sofa::linearalgebra::CompressedRowSparseMatrix<sofa::type::Mat<3, 3, SReal> >::Name();
msg_advice(header) << "Template is empty\n"
<< "By default " << this->getClassName() << " uses blocks with a single scalar (to handle all cases of simulations).\n"
<< "If you are using only 3D DOFs, you may consider using blocks of Matrix3 to speedup the calculations.\n"
<< "If it is the case, add template=\"" << blocksType << "\" to this object in your scene\n"
<< "Otherwise, if you want to disable this message, add " << "template=\"" << this->getTemplateName() << "\" " << ".";
}
if (arg->getAttribute("savingMatrixToFile"))
{
msg_warning() << "It is no longer possible to export the linear system matrix from within " << this->getClassName() << ". Instead, use the component GlobalSystemMatrixExporter (from the SofaMatrix plugin).";
}
}
template<class TMatrix, class TVector, class TThreadManager>
void SparseLDLSolver<TMatrix,TVector,TThreadManager>::solve (Matrix& M, Vector& z, Vector& r)
{
SCOPED_TIMER_VARNAME(solveTimer, "solve");
Inherit::solve_cpu(z.ptr(), r.ptr(), (InvertData *) this->getMatrixInvertData(&M));
}
template <class TMatrix, class TVector, class TThreadManager>
bool SparseLDLSolver<TMatrix, TVector, TThreadManager>::factorize(
Matrix& M, InvertData * invertData)
{
Mfiltered.copyNonZeros(M);
Mfiltered.compress();
int n = M.colSize();
if (n == 0)
{
showInvalidSystemMessage("null size");
return true;
}
int * M_colptr = (int *)Mfiltered.getRowBegin().data();
int * M_rowind = (int *)Mfiltered.getColsIndex().data();
Real * M_values = (Real *)Mfiltered.getColsValue().data();
if (M_colptr == nullptr || M_rowind == nullptr || M_values == nullptr)
{
showInvalidSystemMessage("invalid matrix data structure");
return true;
}
if (Mfiltered.getRowBegin().size() < (size_t)n)
{
showInvalidSystemMessage("size mismatch");
return true;
}
Inherit::factorize(n,M_colptr,M_rowind,M_values, invertData);
numStep++;
return false;
}
template <class TMatrix, class TVector, class TThreadManager>
void SparseLDLSolver<TMatrix, TVector, TThreadManager>::showInvalidSystemMessage(const std::string& reason) const
{
msg_warning() << "Invalid Linear System to solve (" << reason << "). Please insure that there is enough constraints (not rank deficient).";
}
template<class TMatrix, class TVector, class TThreadManager>
void SparseLDLSolver<TMatrix,TVector,TThreadManager>::invert(Matrix& M)
{
factorize(M, (InvertData *) this->getMatrixInvertData(&M));
}
template <class TMatrix, class TVector, class TThreadManager>
bool SparseLDLSolver<TMatrix, TVector, TThreadManager>::doAddJMInvJtLocal(ResMatrixType* result, const JMatrixType* J, SReal fact, InvertData* data)
{
if (!this->isComponentStateValid())
{
return true;
}
/*
J * M^-1 * J^T = J * (L*D*L^T)^-1 * J^t
= (J * (L^T)^-1) * D^-1 * (L^-1 * J^T)
= (L^-1 * J^T)^T * D^-1 * (L^-1 * J^T)
*/
if (J->rowSize() == 0)
{
return true;
}
Jlocal2global.clear();
Jlocal2global.reserve(J->rowSize());
for (auto jit = J->begin(), jitend = J->end(); jit != jitend; ++jit)
{
sofa::SignedIndex l = jit->first;
Jlocal2global.push_back(l);
}
if (Jlocal2global.empty())
{
return true;
}
const unsigned int JlocalRowSize = (unsigned int)Jlocal2global.size();
const simulation::ForEachExecutionPolicy execution = this->d_parallelInverseProduct.getValue() ?
simulation::ForEachExecutionPolicy::PARALLEL :
simulation::ForEachExecutionPolicy::SEQUENTIAL;
simulation::TaskScheduler* taskScheduler = simulation::MainTaskSchedulerFactory::createInRegistry();
assert(taskScheduler);
JLinv.clear();
JLinv.resize(J->rowSize(), data->n);
JLinvDinv.resize(J->rowSize(), data->n);
// copy J in to JLinv taking into account the permutation
unsigned int localRow = 0;
for (auto jit = J->begin(), jitend = J->end(); jit != jitend; ++jit, ++localRow)
{
Real* line = JLinv[localRow];
for (auto it = jit->second.begin(), i2end = jit->second.end(); it != i2end; ++it)
{
int col = data->invperm[it->first];
Real val = it->second;
line[col] = val;
}
}
{
SCOPED_TIMER("LowerSystem");
simulation::forEachRange(execution, *taskScheduler, 0u, JlocalRowSize,
[&data, this](const auto& range)
{
SCOPED_TIMER("Lower");
for (auto i = range.start; i != range.end; ++i)
{
Real* line = JLinv[i];
sofa::linearalgebra::solveLowerUnitriangularSystemCSR(data->n, line, line, data->LT_colptr.data(), data->LT_rowind.data(), data->LT_values.data());
}
});
}
{
SCOPED_TIMER("Diagonal");
simulation::forEachRange(execution, *taskScheduler, 0u, JlocalRowSize,
[&data, this](const auto& range)
{
for (auto i = range.start; i != range.end; ++i)
{
Real* lineD = JLinv[i];
Real* lineM = JLinvDinv[i];
sofa::linearalgebra::solveDiagonalSystemUsingInvertedValues(data->n, lineD, lineM, data->invD.data());
}
});
}
const auto nbTriplets = JlocalRowSize * (JlocalRowSize+1) / 2;
std::vector<Triplet> tripletsBuffer(nbTriplets);
SCOPED_TIMER("UpperSystem");
std::mutex mutex;
// Distribution of the tasks according to the number of triplets, i.e. the
// number of elements in a triangular matrix
simulation::forEachRange(execution, *taskScheduler, 0u, nbTriplets,
[&data, this, fact, &mutex, result, &tripletsBuffer](const auto& range)
{
{
SCOPED_TIMER("UpperRange");
for (auto r = range.start; r != range.end; ++r)
{
//convert a triangular matrix (flat) index to row and column coordinates
sofa::Index i, j;
linearalgebra::computeRowColumnCoordinateFromIndexInLowerTriangularMatrix(r, i, j);
auto& [row, col, value] = tripletsBuffer[r];
row = Jlocal2global[j];
col = Jlocal2global[i];
Real* lineI = JLinv[i];
Real* lineJ = JLinvDinv[j];
value = 0;
for (int k = 0; k < data->n; ++k)
{
value += lineJ[k] * lineI[k];
}
value *= fact;
}
}
std::lock_guard guard(mutex);
SCOPED_TIMER("Assembling");
for (auto r = range.start; r != range.end; ++r)
{
const auto& [row, col, value] = tripletsBuffer[r];
result->add(row, col, value);
if (row != col)
{
result->add(col, row, value);
}
}
});
return true;
}
// Default implementation of Multiply the inverse of the system matrix by the transpose of the given matrix, and multiply the result with the given matrix J
template<class TMatrix, class TVector, class TThreadManager>
bool SparseLDLSolver<TMatrix,TVector,TThreadManager>::addJMInvJtLocal(TMatrix * M, ResMatrixType * result,const JMatrixType * J, SReal fact)
{
InvertData* data = (InvertData*)this->getMatrixInvertData(M);
return doAddJMInvJtLocal(result, J, fact, data);
}
} // namespace sofa::component::linearsolver::direct