5#if HAVE_SUITESPARSE_CHOLMOD
12#include <dune/istl/solverfactory.hh>
13#include <dune/istl/foreach.hh>
34 const NoIgnore& operator[](std::size_t)
const {
return *
this; }
35 explicit operator bool()
const {
return false; }
36 static constexpr std::size_t size() {
return 0; }
41 template<
class BlockedVector,
class FlatVector>
42 void copyToFlatVector(
const BlockedVector& blockedVector, FlatVector& flatVector)
46 flatVector.resize(len);
49 flatVector[offset] = entry;
54 template<
class FlatVector>
55 void copyToFlatVector(
const NoIgnore&, FlatVector&)
61 template<
class FlatVector,
class BlockedVector>
62 void copyToBlockedVector(
const FlatVector& flatVector, BlockedVector& blockedVector)
65 entry = flatVector[offset];
77class Cholmod :
public InverseOperator<Vector, Vector>
99 cholmod_free_factor(&L_, &c_);
104 Cholmod(
const Cholmod&) =
delete;
105 Cholmod& operator=(
const Cholmod&) =
delete;
110 void apply (Vector& x, Vector& b, [[maybe_unused]]
double reduction, InverseOperatorResult& res)
120 void apply(Vector& x, Vector& b, InverseOperatorResult& res)
128 if (x.size() != b.size())
129 DUNE_THROW(Exception,
"Error in apply(): sizes of x and b do not match!");
132 auto b2 = std::make_unique<double[]>(L_->n);
133 auto x2 = std::make_unique<double[]>(L_->n);
139 if ( subIndices_.empty() )
140 bp[ flatIndex ] = entry;
143 bp[ subIndices_[ flatIndex ] ] = entry;
147 auto b3 = make_cholmod_dense(cholmod_allocate_dense(L_->n, 1, L_->n, CHOLMOD_REAL, &c_), &c_);
149 auto b4 =
static_cast<double*
>(b3->x);
150 std::copy(b2.get(), b2.get() + L_->n, b4);
153 auto x3 = make_cholmod_dense(cholmod_solve(CHOLMOD_A, L_, b3.get(), &c_), &c_);
155 auto xp =
static_cast<double*
>(x3->x);
159 if ( subIndices_.empty() )
160 entry = xp[ flatIndex ];
163 entry = xp[ subIndices_[ flatIndex ] ];
168 res.converged =
true;
177 template<
class Matrix>
178 void setMatrix(
const Matrix& matrix)
180 const Impl::NoIgnore* noIgnore =
nullptr;
181 setMatrix(matrix, noIgnore);
198 template<
class Matrix,
class Ignore>
199 void setMatrix(
const Matrix& matrix,
const Ignore* ignore)
203 int numberOfIgnoredDofs = 0;
206 auto [flatRows,flatCols] =
flatMatrixForEach( matrix, [&](
auto&& ,
auto&& flatRowIndex,
auto&& flatColIndex){
207 if( flatRowIndex <= flatColIndex )
211 std::vector<bool> flatIgnore;
215 Impl::copyToFlatVector(*ignore,flatIgnore);
216 numberOfIgnoredDofs = std::count(flatIgnore.begin(),flatIgnore.end(),
true);
220 int N = flatRows - numberOfIgnoredDofs;
234 const auto deleter = [c = &this->c_](
auto* p) {
235 cholmod_free_sparse(&p, c);
237 auto M = std::unique_ptr<cholmod_sparse, decltype(deleter)>(
238 cholmod_allocate_sparse(N,
249 int* Ap =
static_cast<int*
>(M->p);
250 int* Ai =
static_cast<int*
>(M->i);
251 double* Ax =
static_cast<double*
>(M->x);
259 std::size_t subIndexCounter = 0;
261 for ( std::size_t i=0; i<flatRows; i++ )
263 if ( not flatIgnore[ i ] )
265 subIndices_[ i ] = subIndexCounter++;
272 flatMatrixForEach(matrix, [&](
auto&& ,
auto&& flatRowIndex,
auto&& flatColIndex){
275 if ( ignore and ( flatIgnore[flatRowIndex] or flatIgnore[flatColIndex] ) )
279 if ( flatRowIndex > flatColIndex )
283 auto idx = ignore ? subIndices_[flatRowIndex] : flatRowIndex;
290 for (
int i=0; i<N; i++ )
296 std::vector<std::size_t> rowPosition(N,0);
299 flatMatrixForEach(matrix, [&](
auto&& entry,
auto&& flatRowIndex,
auto&& flatColIndex){
302 if ( ignore and ( flatIgnore[flatRowIndex] or flatIgnore[flatColIndex] ) )
306 if ( flatRowIndex > flatColIndex )
310 auto rowIdx = ignore ? subIndices_[flatRowIndex] : flatRowIndex;
311 auto colIdx = ignore ? subIndices_[flatColIndex] : flatColIndex;
312 auto rowStart = Ap[rowIdx];
313 auto rowPos = rowPosition[rowIdx];
314 Ai[ rowStart + rowPos ] = colIdx;
315 Ax[ rowStart + rowPos ] = entry;
316 rowPosition[rowIdx]++;
321 L_ = cholmod_analyze(M.get(), &c_);
324 cholmod_factorize(M.get(), L_, &c_);
329 return SolverCategory::Category::sequential;
337 cholmod_common& cholmodCommonObject()
347 cholmod_factor& cholmodFactor()
357 const cholmod_factor& cholmodFactor()
const
365 auto make_cholmod_dense(cholmod_dense* x, cholmod_common* c)
367 const auto deleter = [c](
auto* p) {
368 cholmod_free_dense(&p, c);
370 return std::unique_ptr<cholmod_dense, decltype(deleter)>(x, deleter);
374 cholmod_factor* L_ =
nullptr;
377 bool nIsZero_ =
false;
380 std::vector<std::size_t> subIndices_;
383 struct CholmodCreator{
384 template<
class F>
struct isValidBlock : std::false_type{};
385 template<
int k>
struct isValidBlock<FieldVector<double,k>> : std::true_type{};
386 template<
int k>
struct isValidBlock<FieldVector<float,k>> : std::true_type{};
388 template<
class TL,
typename M>
389 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
390 typename Dune::TypeListElement<2, TL>::type>>
392 std::enable_if_t<isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
394 using D =
typename Dune::TypeListElement<1, TL>::type;
395 auto solver = std::make_shared<Dune::Cholmod<D>>();
396 solver->setMatrix(mat);
401 template<
typename TL,
typename M>
402 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
403 typename Dune::TypeListElement<2, TL>::type>>
405 std::enable_if_t<!isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
407 DUNE_THROW(UnsupportedType,
"Unsupported Type in Cholmod");
410 DUNE_REGISTER_DIRECT_SOLVER(
"cholmod", Dune::CholmodCreator());
Implementation of the BCRSMatrix class.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Hierarchical structure of string parameters.
Definition: parametertree.hh:37
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
Dune namespace.
Definition: alignedallocator.hh:13
std::pair< std::size_t, std::size_t > flatMatrixForEach(Matrix &&matrix, F &&f, std::size_t rowOffset=0, std::size_t colOffset=0)
Traverse a blocked matrix and call a functor at each scalar entry.
Definition: foreach.hh:132
std::size_t flatVectorForEach(Vector &&vector, F &&f, std::size_t offset=0)
Traverse a blocked vector and call a functor at each scalar entry.
Definition: foreach.hh:95
Define general, extensible interface for inverse operators.
Category
Definition: solvercategory.hh:23