5#ifndef DUNE_ISTL_SUPERMATRIX_HH
6#define DUNE_ISTL_SUPERMATRIX_HH
17#include <dune/istl/bccsmatrixinitializer.hh>
19#include "superlufunctions.hh"
25 struct SuperMatrixCreateSparseChooser
29 struct SuperMatrixPrinter
32#if __has_include("slu_sdefs.h")
34 struct SuperMatrixCreateSparseChooser<float>
36 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
37 float *values,
int *rowindex,
int* colindex,
38 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
40 sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
46 struct SuperMatrixPrinter<float>
48 static void print(
char* name, SuperMatrix* mat)
50 sPrint_CompCol_Matrix(name, mat);
55#if __has_include("slu_ddefs.h")
57 struct SuperMatrixCreateSparseChooser<double>
59 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
60 double *values,
int *rowindex,
int* colindex,
61 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
63 dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
69 struct SuperMatrixPrinter<double>
71 static void print(
char* name, SuperMatrix* mat)
73 dPrint_CompCol_Matrix(name, mat);
78#if __has_include("slu_cdefs.h")
80 struct SuperMatrixCreateSparseChooser<
std::complex<float> >
82 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
83 std::complex<float> *values,
int *rowindex,
int* colindex,
84 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
86 cCreate_CompCol_Matrix(mat, n, m, offset,
reinterpret_cast< ::complex*
>(values),
87 rowindex, colindex, stype, dtype, mtype);
92 struct SuperMatrixPrinter<
std::complex<float> >
94 static void print(
char* name, SuperMatrix* mat)
96 cPrint_CompCol_Matrix(name, mat);
101#if __has_include("slu_zdefs.h")
103 struct SuperMatrixCreateSparseChooser<
std::complex<double> >
105 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
106 std::complex<double> *values,
int *rowindex,
int* colindex,
107 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
109 zCreate_CompCol_Matrix(mat, n, m, offset,
reinterpret_cast<doublecomplex*
>(values),
110 rowindex, colindex, stype, dtype, mtype);
115 struct SuperMatrixPrinter<
std::complex<double> >
117 static void print(
char* name, SuperMatrix* mat)
119 zPrint_CompCol_Matrix(name, mat);
125 struct BaseGetSuperLUType
127 static const Dtype_t type;
131 struct GetSuperLUType
135 const Dtype_t BaseGetSuperLUType<T>::type =
136 std::is_same<T,float>::value ? SLU_S :
137 ( std::is_same<T,std::complex<double> >::value ? SLU_Z :
138 ( std::is_same<T,std::complex<float> >::value ? SLU_C : SLU_D ));
141 struct GetSuperLUType<double>
142 :
public BaseGetSuperLUType<double>
144 typedef double float_type;
148 struct GetSuperLUType<float>
149 :
public BaseGetSuperLUType<float>
151 typedef float float_type;
155 struct GetSuperLUType<
std::complex<double> >
156 :
public BaseGetSuperLUType<std::complex<double> >
158 typedef double float_type;
162 struct GetSuperLUType<
std::complex<float> >
163 :
public BaseGetSuperLUType<std::complex<float> >
165 typedef float float_type;
178 struct SuperMatrixInitializer
184 template<
class M,
class X,
class TM,
class TD,
class T1>
185 class SeqOverlappingSchwarz;
187 template<
class T,
bool flag>
188 struct SeqOverlappingSchwarzAssemblerHelper;
193 template<
class B,
class TA>
195 :
public ISTL::Impl::BCCSMatrix<typename BCRSMatrix<B,TA>::field_type, int>
197 template<
class M,
class X,
class TM,
class TD,
class T1>
199 friend struct SuperMatrixInitializer<
BCRSMatrix<B,TA> >;
204 friend struct SeqOverlappingSchwarzAssemblerHelper<
SuperLU<
Matrix>, true>;
221 if (this->N_+this->M_*this->Nnz_ != 0)
226 operator SuperMatrix&()
232 operator const SuperMatrix&()
const
239 if (this->N_ + this->M_ + this->Nnz_!=0)
243 this->N_ = MatrixDimension<Matrix>::rowdim(mat);
244 this->M_ = MatrixDimension<Matrix>::coldim(mat);
245 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(*
this);
247 copyToBCCSMatrix(initializer, mat);
249 SuperMatrixCreateSparseChooser<typename Matrix::field_type>
250 ::create(&A, this->N_, this->M_, this->colstart[this->N_],
251 this->values,this->rowindex, this->colstart, SLU_NC,
252 static_cast<Dtype_t
>(GetSuperLUType<typename Matrix::field_type>::type), SLU_GE);
256 SuperLUMatrix<BCRSMatrix<B,TA> >& operator=(
const SuperLUMatrix <BCRSMatrix<B,TA> >& mat)
258 if (this->N_ + this->M_ + this->Nnz_!=0)
261 using Matrix = BCRSMatrix<B,TA>;
262 this->N_ = MatrixDimension<Matrix>::rowdim(mat);
263 this->M_ = MatrixDimension<Matrix>::coldim(mat);
264 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(*
this);
266 copyToBCCSMatrix(initializer, mat);
268 SuperMatrixCreateSparseChooser<B>
269 ::create(&A, this->N_, this->M_, this->colstart[this->N_],
270 this->values,this->rowindex, this->colstart, SLU_NC,
271 static_cast<Dtype_t
>(GetSuperLUType<B>::type), SLU_GE);
283 if(this->N_+this->M_+this->Nnz_!=0)
285 this->N_=mrs.size()*MatrixDimension<typename Matrix::block_type>::rowdim(*(mat[0].begin()));
286 this->M_=mrs.size()*MatrixDimension<typename Matrix::block_type>::coldim(*(mat[0].begin()));
287 SuperMatrixInitializer<Matrix> initializer(*
this);
289 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<
Matrix,std::set<std::size_t> >(mat,mrs));
295 this->N_=MatrixDimension<Matrix>::rowdim(mat);
296 this->M_=MatrixDimension<Matrix>::coldim(mat);
297 SuperMatrixInitializer<Matrix> initializer(*
this);
299 copyToBCCSMatrix(initializer, mat);
305 ISTL::Impl::BCCSMatrix<typename BCRSMatrix<B,TA>::field_type,
int>::free();
306 SUPERLU_FREE(A.Store);
312 template<
class B,
class A>
313 class SuperMatrixInitializer<BCRSMatrix<B,A> >
314 :
public ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<B,A>, int>
316 template<
class I,
class S,
class D>
317 friend class OverlappingSchwarzInitializer;
319 typedef BCRSMatrix<B,A> Matrix;
322 SuperMatrixInitializer(SuperLUMatrix& lum) : ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<B,A>, int>(lum)
326 SuperMatrixInitializer() : ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<B,A>, int>()
329 void createMatrix()
const override
331 ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<B,A>,
int>::createMatrix();
332 SuperMatrixCreateSparseChooser<typename Matrix::field_type>
333 ::create(&slumat->A, slumat->N_, slumat->M_, slumat->colstart[this->cols],
334 slumat->values,slumat->rowindex, slumat->colstart, SLU_NC,
335 static_cast<Dtype_t
>(GetSuperLUType<typename Matrix::field_type>::type), SLU_GE);
338 SuperLUMatrix* slumat;
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:466
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:497
A generic dynamic dense matrix.
Definition: matrix.hh:561
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:755
SuperLUMatrix(const Matrix &mat)
Constructor that initializes the data.
Definition: supermatrix.hh:212
virtual void setMatrix(const Matrix &mat)
Initialize data from given matrix.
Definition: supermatrix.hh:293
BCRSMatrix< B, TA > Matrix
The type of the matrix to convert.
Definition: supermatrix.hh:202
virtual void setMatrix(const Matrix &mat, const std::set< std::size_t > &mrs)
Initialize data from a given set of matrix rows and columns.
Definition: supermatrix.hh:281
virtual ~SuperLUMatrix()
Destructor.
Definition: supermatrix.hh:219
void free() override
free allocated space.
Definition: supermatrix.hh:303
SuperLu Solver.
Definition: superlu.hh:271
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.
Dune namespace.
Definition: alignedallocator.hh:13
Utility class for converting an ISTL Matrix into a SuperLU Matrix.
Definition: supermatrix.hh:175
Traits for type conversions and type information.