3#if HAVE_SUITESPARSE_CHOLMOD
10#include <dune/istl/solverfactory.hh>
11#include <dune/istl/foreach.hh>
32 const NoIgnore& operator[](std::size_t)
const {
return *
this; }
33 explicit operator bool()
const {
return false; }
34 static constexpr std::size_t size() {
return 0; }
39 template<
class BlockedVector,
class FlatVector>
40 void copyToFlatVector(
const BlockedVector& blockedVector, FlatVector& flatVector)
44 flatVector.resize(len);
47 flatVector[offset] = entry;
52 template<
class FlatVector>
53 void copyToFlatVector(
const NoIgnore&, FlatVector&)
59 template<
class FlatVector,
class BlockedVector>
60 void copyToBlockedVector(
const FlatVector& flatVector, BlockedVector& blockedVector)
63 entry = flatVector[offset];
75class Cholmod :
public InverseOperator<Vector, Vector>
97 cholmod_free_factor(&L_, &c_);
102 Cholmod(
const Cholmod&) =
delete;
103 Cholmod& operator=(
const Cholmod&) =
delete;
108 void apply (Vector& x, Vector& b, [[maybe_unused]]
double reduction, InverseOperatorResult& res)
118 void apply(Vector& x, Vector& b, InverseOperatorResult& res)
126 if (x.size() != b.size())
127 DUNE_THROW(Exception,
"Error in apply(): sizes of x and b do not match!");
130 auto b2 = std::make_unique<double[]>(L_->n);
131 auto x2 = std::make_unique<double[]>(L_->n);
137 if ( subIndices_.empty() )
138 bp[ flatIndex ] = entry;
141 bp[ subIndices_[ flatIndex ] ] = entry;
145 auto b3 = make_cholmod_dense(cholmod_allocate_dense(L_->n, 1, L_->n, CHOLMOD_REAL, &c_), &c_);
147 auto b4 =
static_cast<double*
>(b3->x);
148 std::copy(b2.get(), b2.get() + L_->n, b4);
151 auto x3 = make_cholmod_dense(cholmod_solve(CHOLMOD_A, L_, b3.get(), &c_), &c_);
153 auto xp =
static_cast<double*
>(x3->x);
157 if ( subIndices_.empty() )
158 entry = xp[ flatIndex ];
161 entry = xp[ subIndices_[ flatIndex ] ];
166 res.converged =
true;
175 template<
class Matrix>
176 void setMatrix(
const Matrix& matrix)
178 const Impl::NoIgnore* noIgnore =
nullptr;
179 setMatrix(matrix, noIgnore);
196 template<
class Matrix,
class Ignore>
197 void setMatrix(
const Matrix& matrix,
const Ignore* ignore)
201 int numberOfIgnoredDofs = 0;
204 auto [flatRows,flatCols] =
flatMatrixForEach( matrix, [&](
auto&& ,
auto&& flatRowIndex,
auto&& flatColIndex){
205 if( flatRowIndex <= flatColIndex )
209 std::vector<bool> flatIgnore;
213 Impl::copyToFlatVector(*ignore,flatIgnore);
214 numberOfIgnoredDofs = std::count(flatIgnore.begin(),flatIgnore.end(),
true);
218 int N = flatRows - numberOfIgnoredDofs;
232 const auto deleter = [c = &this->c_](
auto* p) {
233 cholmod_free_sparse(&p, c);
235 auto M = std::unique_ptr<cholmod_sparse, decltype(deleter)>(
236 cholmod_allocate_sparse(N,
247 int* Ap =
static_cast<int*
>(M->p);
248 int* Ai =
static_cast<int*
>(M->i);
249 double* Ax =
static_cast<double*
>(M->x);
257 std::size_t subIndexCounter = 0;
259 for ( std::size_t i=0; i<flatRows; i++ )
261 if ( not flatIgnore[ i ] )
263 subIndices_[ i ] = subIndexCounter++;
270 flatMatrixForEach(matrix, [&](
auto&& ,
auto&& flatRowIndex,
auto&& flatColIndex){
273 if ( ignore and ( flatIgnore[flatRowIndex] or flatIgnore[flatColIndex] ) )
277 if ( flatRowIndex > flatColIndex )
281 auto idx = ignore ? subIndices_[flatRowIndex] : flatRowIndex;
288 for (
int i=0; i<N; i++ )
294 std::vector<std::size_t> rowPosition(N,0);
297 flatMatrixForEach(matrix, [&](
auto&& entry,
auto&& flatRowIndex,
auto&& flatColIndex){
300 if ( ignore and ( flatIgnore[flatRowIndex] or flatIgnore[flatColIndex] ) )
304 if ( flatRowIndex > flatColIndex )
308 auto rowIdx = ignore ? subIndices_[flatRowIndex] : flatRowIndex;
309 auto colIdx = ignore ? subIndices_[flatColIndex] : flatColIndex;
310 auto rowStart = Ap[rowIdx];
311 auto rowPos = rowPosition[rowIdx];
312 Ai[ rowStart + rowPos ] = colIdx;
313 Ax[ rowStart + rowPos ] = entry;
314 rowPosition[rowIdx]++;
319 L_ = cholmod_analyze(M.get(), &c_);
322 cholmod_factorize(M.get(), L_, &c_);
327 return SolverCategory::Category::sequential;
335 cholmod_common& cholmodCommonObject()
343 auto make_cholmod_dense(cholmod_dense* x, cholmod_common* c)
345 const auto deleter = [c](
auto* p) {
346 cholmod_free_dense(&p, c);
348 return std::unique_ptr<cholmod_dense, decltype(deleter)>(x, deleter);
352 cholmod_factor* L_ =
nullptr;
355 bool nIsZero_ =
false;
358 std::vector<std::size_t> subIndices_;
361 struct CholmodCreator{
362 template<
class F>
struct isValidBlock : std::false_type{};
363 template<
int k>
struct isValidBlock<FieldVector<double,k>> : std::true_type{};
364 template<
int k>
struct isValidBlock<FieldVector<float,k>> : std::true_type{};
366 template<
class TL,
typename M>
367 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
368 typename Dune::TypeListElement<2, TL>::type>>
370 std::enable_if_t<isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
372 using D =
typename Dune::TypeListElement<1, TL>::type;
373 auto solver = std::make_shared<Dune::Cholmod<D>>();
374 solver->setMatrix(mat);
379 template<
typename TL,
typename M>
380 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
381 typename Dune::TypeListElement<2, TL>::type>>
383 std::enable_if_t<!isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
385 DUNE_THROW(UnsupportedType,
"Unsupported Type in Cholmod");
388 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:35
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:216
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
Dune namespace.
Definition: alignedallocator.hh:11
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:129
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:92
Define general, extensible interface for inverse operators.
Category
Definition: solvercategory.hh:21