Dune Core Modules (2.7.0)

cholmod.hh
1#pragma once
2
8#include <dune/istl/solverfactory.hh>
9
10#include <vector>
11#include <memory>
12
13#include <cholmod.h>
14
15namespace Dune {
16
17namespace Impl{
18
27 struct NoIgnore
28 {
29 const NoIgnore& operator[](std::size_t) const { return *this; }
30 explicit operator bool() const { return false; }
31 std::size_t count() const { return 0; }
32 };
33
34} //namespace Impl
35
36template<class T>
37class Cholmod;
38
43template<class T, class A, int k>
44class Cholmod<BlockVector<FieldVector<T,k>, A>>
45 : public InverseOperator<BlockVector<FieldVector<T,k>, A>, BlockVector<FieldVector<T,k>, A>>
46{
47public:
48
49 // type of unknown
51 // type of rhs
53
54
61 {
62 cholmod_start(&c_);
63 }
64
71 {
72 if (L_)
73 cholmod_free_factor(&L_, &c_);
74 cholmod_finish(&c_);
75 }
76
77 // forbid copying to avoid freeing memory twice
78 Cholmod(const Cholmod&) = delete;
79 Cholmod& operator=(const Cholmod&) = delete;
80
81
84 void apply (X& x, B& b, double reduction, InverseOperatorResult& res)
85 {
86 DUNE_UNUSED_PARAMETER(reduction);
87 apply(x,b,res);
88 }
89
95 void apply(X& x, B& b, InverseOperatorResult& res)
96 {
97 // do nothing if N=0
98 if ( nIsZero_ )
99 {
100 return;
101 }
102
103 if (x.size() != b.size())
104 DUNE_THROW(Exception, "Error in apply(): sizes of x and b do not match!");
105
106 const auto& blocksize = k;
107
108 // cast to double array
109 auto b2 = std::make_unique<double[]>(L_->n);
110 auto x2 = std::make_unique<double[]>(L_->n);
111
112 // copy to cholmod
113 auto bp = b2.get();
114 if (inverseSubIndices_.empty()) // no ignore field given
115 {
116 // simply copy all values
117 for (const auto& bi : b)
118 for (const auto& bii : bi)
119 *bp++ = bii;
120 }
121 else // use the mapping from not ignored entries
122 {
123 // iterate over inverseSubIndices and resolve the block indices
124 for (const auto& idx : inverseSubIndices_)
125 *bp++ = b[ idx / blocksize ][ idx % blocksize ];
126 }
127
128 // create a cholmod dense object
129 auto b3 = make_cholmod_dense(cholmod_allocate_dense(L_->n, 1, L_->n, CHOLMOD_REAL, &c_), &c_);
130 // cast because void-ptr
131 auto b4 = static_cast<double*>(b3->x);
132 std::copy(b2.get(), b2.get() + L_->n, b4);
133
134 // solve for a cholmod x object
135 auto x3 = make_cholmod_dense(cholmod_solve(CHOLMOD_A, L_, b3.get(), &c_), &c_);
136 // cast because void-ptr
137 auto xp = static_cast<double*>(x3->x);
138
139 // copy into x
140 if (inverseSubIndices_.empty()) // no ignore field given
141 {
142 // simply copy all values
143 for (int i=0, s=x.size(); i<s; i++)
144 for (int ii=0, ss=x[i].size(); ii<ss; ii++)
145 x[i][ii] = *xp++;
146 }
147 else // use the mapping from not ignored entries
148 {
149 // iterate over inverseSubIndices and resolve the block indices
150 for (const auto& idx : inverseSubIndices_)
151 x[ idx / blocksize ][ idx % blocksize ] = *xp++;
152 }
153
154 // statistics for a direct solver
155 res.iterations = 1;
156 res.converged = true;
157 }
158
159
166 {
167 const Impl::NoIgnore* noIgnore = nullptr;
168 setMatrix(matrix, noIgnore);
169 }
170
185 template<class Ignore>
186 void setMatrix(const BCRSMatrix<FieldMatrix<T,k,k>>& matrix, const Ignore* ignore)
187 {
188
189 const auto blocksize = k;
190
191 // Total number of rows
192 int N = blocksize * matrix.N();
193 if ( ignore )
194 N -= ignore->count();
195
196 nIsZero_ = (N <= 0);
197
198 if ( nIsZero_ )
199 {
200 return;
201 }
202
203 // number of nonzeroes
204 const int nnz = blocksize * blocksize * matrix.nonzeroes();
205 // number of diagonal entries
206 const int nDiag = blocksize * blocksize * matrix.N();
207 // number of nonzeroes in the dialgonal
208 const int nnzDiag = (blocksize * (blocksize+1)) / 2 * matrix.N();
209 // number of nonzeroes in triangular submatrix (including diagonal)
210 const int nnzTri = (nnz - nDiag) / 2 + nnzDiag;
211
212 /*
213 * CHOLMOD uses compressed-column sparse matrices, but for symmetric
214 * matrices this is the same as the compressed-row sparse matrix used
215 * by DUNE. So we can just store Mᵀ instead of M (as M = Mᵀ).
216 */
217 const auto deleter = [c = &this->c_](auto* p) {
218 cholmod_free_sparse(&p, c);
219 };
220 auto M = std::unique_ptr<cholmod_sparse, decltype(deleter)>(
221 cholmod_allocate_sparse(N, // # rows
222 N, // # cols
223 nnzTri, // # of nonzeroes
224 1, // indices are sorted ( 1 = true)
225 1, // matrix is "packed" ( 1 = true)
226 -1, // stype of matrix ( -1 = cosider the lower part only )
227 CHOLMOD_REAL, // xtype of matrix ( CHOLMOD_REAL = single array, no complex numbers)
228 &c_ // cholmod_common ptr
229 ), deleter);
230
231 // copy the data of BCRS matrix to Cholmod Sparse matrix
232 int* Ap = static_cast<int*>(M->p);
233 int* Ai = static_cast<int*>(M->i);
234 double* Ax = static_cast<double*>(M->x);
235
236 // vector mapping all indices in flat order to the not ignored indices
237 std::vector<std::size_t> subIndices;
238
239 if ( ignore )
240 {
241 // init the mappings
242 inverseSubIndices_.resize(N); // size = number of not ignored entries
243 subIndices.resize(matrix.M()*blocksize); // size = number of all entries
244
245 std::size_t j=0;
246 for( std::size_t block=0; block<matrix.N(); block++ )
247 {
248 for( std::size_t i=0; i<blocksize; i++ )
249 {
250 if( not (*ignore)[block][i] )
251 {
252 subIndices[ block*blocksize + i ] = j;
253 inverseSubIndices_[j] = block*blocksize + i;
254 j++;
255 }
256 }
257 }
258 }
259
260 // Copy the data to the CHOLMOD array
261 int n = 0;
262 for (auto rowIt = matrix.begin(); rowIt != matrix.end(); rowIt++)
263 {
264 const auto row = rowIt.index();
265 for (std::size_t i=0; i<blocksize; i++)
266 {
267 if( ignore and (*ignore)[row][i] )
268 continue;
269
270 // col start
271 *Ap++ = n;
272
273 for (auto colIt = rowIt->begin(); colIt != rowIt->end(); ++colIt)
274 {
275 const auto col = colIt.index();
276
277 // are we already in the lower part?
278 if (col < row)
279 continue;
280
281 for (auto j = (row == col ? i : 0); j<blocksize; j++)
282 {
283 if( ignore and (*ignore)[col][j] )
284 continue;
285
286 const auto jj = ignore ? subIndices[ blocksize*col + j ] : blocksize*col + j;
287
288 // set the current index and entry
289 *Ai++ = jj;
290 *Ax++ = (*colIt)[i][j];
291 // increase number of set values
292 n++;
293 }
294 }
295 }
296 }
297 // set last col start
298 *Ap = n;
299
300 // Now analyse the pattern and optimal row order
301 L_ = cholmod_analyze(M.get(), &c_);
302
303 // Do the factorization (this may take some time)
304 cholmod_factorize(M.get(), L_, &c_);
305 }
306
308 {
309 return SolverCategory::Category::sequential;
310 }
311
312private:
313
314 // create a destrucable unique_ptr
315 auto make_cholmod_dense(cholmod_dense* x, cholmod_common* c)
316 {
317 const auto deleter = [c](auto* p) {
318 cholmod_free_dense(&p, c);
319 };
320 return std::unique_ptr<cholmod_dense, decltype(deleter)>(x, deleter);
321 }
322
323 cholmod_common c_;
324 cholmod_factor* L_ = nullptr;
325
326 // indicator for a 0x0 problem (due to ignore dof's)
327 bool nIsZero_ = false;
328
329 // mapping from the not ignored indices in flat order to all indices in flat order
330 // it also holds the info about ignore nodes: if it is empty => no ignore field
331 std::vector<std::size_t> inverseSubIndices_;
332
333};
334
335 struct CholmodCreator{
336 template<class F> struct isValidBlock : std::false_type{};
337 template<int k> struct isValidBlock<FieldVector<double,k>> : std::true_type{};
338 template<int k> struct isValidBlock<FieldVector<float,k>> : std::true_type{};
339
340 template<class TL, typename M>
341 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
342 typename Dune::TypeListElement<2, TL>::type>>
343 operator()(TL /*tl*/, const M& mat, const Dune::ParameterTree& /*config*/,
344 std::enable_if_t<isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
345 {
346 using D = typename Dune::TypeListElement<1, TL>::type;
347 auto solver = std::make_shared<Dune::Cholmod<D>>();
348 solver->setMatrix(mat);
349 return solver;
350 };
351
352 // second version with SFINAE to validate the template parameters of Cholmod
353 template<typename TL, typename M>
354 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
355 typename Dune::TypeListElement<2, TL>::type>>
356 operator() (TL /*tl*/, const M& /*mat*/, const Dune::ParameterTree& /*config*/,
357 std::enable_if_t<!isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
358 {
359 DUNE_THROW(UnsupportedType, "Unsupported Type in Cholmod");
360 };
361 };
362 DUNE_REGISTER_DIRECT_SOLVER("cholmod", Dune::CholmodCreator());
363
364} /* namespace Dune */
Implementation of the BCRSMatrix class.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:425
A vector of blocks with memory management.
Definition: bvector.hh:403
void setMatrix(const BCRSMatrix< FieldMatrix< T, k, k > > &matrix, const Ignore *ignore)
Set matrix and ignore nodes.
Definition: cholmod.hh:186
void apply(X &x, B &b, double reduction, InverseOperatorResult &res)
simple forward to apply(X&, Y&, InverseOperatorResult&)
Definition: cholmod.hh:84
void setMatrix(const BCRSMatrix< FieldMatrix< T, k, k > > &matrix)
Set matrix without ignore nodes.
Definition: cholmod.hh:165
~Cholmod()
Destructor.
Definition: cholmod.hh:70
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: cholmod.hh:307
Cholmod()
Default constructor.
Definition: cholmod.hh:60
void apply(X &x, B &b, InverseOperatorResult &res)
solve the linear system Ax=b (possibly with respect to some ignore field)
Definition: cholmod.hh:95
Base class for Dune-Exceptions.
Definition: exceptions.hh:94
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Abstract base class for all solvers.
Definition: solver.hh:97
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.
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:46
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Dune namespace.
Definition: alignedallocator.hh:14
Define general, extensible interface for inverse operators.
Statistics about the application of an inverse operator.
Definition: solver.hh:46
int iterations
Number of iterations.
Definition: solver.hh:65
bool converged
True if convergence criterion has been met.
Definition: solver.hh:71
Category
Definition: solvercategory.hh:21
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)