3#ifndef DUNE_ISTL_SUPERMATRIX_HH
4#define DUNE_ISTL_SUPERMATRIX_HH
15#include"colcompmatrix.hh"
17#include "superlufunctions.hh"
23 struct SuperMatrixCreateSparseChooser
27 struct SuperMatrixPrinter
32 struct SuperMatrixCreateSparseChooser<float>
34 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
35 float *values,
int *rowindex,
int* colindex,
36 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
38 sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
44 struct SuperMatrixPrinter<float>
46 static void print(
char* name, SuperMatrix* mat)
48 sPrint_CompCol_Matrix(name, mat);
55 struct SuperMatrixCreateSparseChooser<double>
57 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
58 double *values,
int *rowindex,
int* colindex,
59 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
61 dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
67 struct SuperMatrixPrinter<double>
69 static void print(
char* name, SuperMatrix* mat)
71 dPrint_CompCol_Matrix(name, mat);
78 struct SuperMatrixCreateSparseChooser<
std::complex<float> >
80 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
81 std::complex<float> *values,
int *rowindex,
int* colindex,
82 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
84 cCreate_CompCol_Matrix(mat, n, m, offset,
reinterpret_cast< ::complex*
>(values),
85 rowindex, colindex, stype, dtype, mtype);
90 struct SuperMatrixPrinter<
std::complex<float> >
92 static void print(
char* name, SuperMatrix* mat)
94 cPrint_CompCol_Matrix(name, mat);
101 struct SuperMatrixCreateSparseChooser<
std::complex<double> >
103 static void create(SuperMatrix *mat,
int n,
int m,
int offset,
104 std::complex<double> *values,
int *rowindex,
int* colindex,
105 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
107 zCreate_CompCol_Matrix(mat, n, m, offset,
reinterpret_cast<doublecomplex*
>(values),
108 rowindex, colindex, stype, dtype, mtype);
113 struct SuperMatrixPrinter<
std::complex<double> >
115 static void print(
char* name, SuperMatrix* mat)
117 zPrint_CompCol_Matrix(name, mat);
123 struct BaseGetSuperLUType
125 static const Dtype_t type;
129 struct GetSuperLUType
133 const Dtype_t BaseGetSuperLUType<T>::type =
134 std::is_same<T,float>::value ? SLU_S :
135 ( std::is_same<T,std::complex<double> >::value ? SLU_Z :
136 ( std::is_same<T,std::complex<float> >::value ? SLU_C : SLU_D ));
139 struct GetSuperLUType<double>
140 :
public BaseGetSuperLUType<double>
142 typedef double float_type;
146 struct GetSuperLUType<float>
147 :
public BaseGetSuperLUType<float>
149 typedef float float_type;
153 struct GetSuperLUType<
std::complex<double> >
154 :
public BaseGetSuperLUType<std::complex<double> >
156 typedef double float_type;
160 struct GetSuperLUType<
std::complex<float> >
161 :
public BaseGetSuperLUType<std::complex<float> >
163 typedef float float_type;
176 struct SuperMatrixInitializer
185 template<
class B,
class TA,
int n,
int m>
189 template<
class M,
class X,
class TM,
class TD,
class T1>
196 friend struct SeqOverlappingSchwarzAssemblerHelper<SuperLU<
Matrix>, true>;
213 if (this->N_+this->M_*this->Nnz_ != 0)
218 operator SuperMatrix&()
224 operator const SuperMatrix&()
const
232 SuperMatrixCreateSparseChooser<B>
233 ::create(&A, this->N_, this->M_, this->colstart[this->N_],
234 this->values,this->rowindex, this->colstart, SLU_NC,
235 static_cast<Dtype_t
>(GetSuperLUType<B>::type), SLU_GE);
239 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& operator=(
const SuperLUMatrix <BCRSMatrix<FieldMatrix<B,n,m>,TA> >& mat)
241 this->ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(mat);
242 SuperMatrixCreateSparseChooser<B>
243 ::create(&A, this->N_, this->M_, this->colstart[this->N_],
244 this->values,this->rowindex, this->colstart, SLU_NC,
245 static_cast<Dtype_t
>(GetSuperLUType<B>::type), SLU_GE);
257 if(this->N_+this->M_+this->Nnz_!=0)
259 this->N_=mrs.size()*n;
260 this->M_=mrs.size()*m;
261 SuperMatrixInitializer<Matrix> initializer(*
this);
271 SuperMatrixInitializer<Matrix> initializer(*
this);
280 SUPERLU_FREE(A.Store);
286 template<
class T,
class A,
int n,
int m>
287 class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
288 :
public ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
290 template<
class I,
class S,
class D>
291 friend class OverlappingSchwarzInitializer;
293 typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
296 SuperMatrixInitializer(SuperLUMatrix& lum) : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >(lum)
300 SuperMatrixInitializer() : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >()
303 virtual void createMatrix()
const
305 ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::createMatrix();
306 SuperMatrixCreateSparseChooser<T>
307 ::create(&slumat->A, slumat->N_, slumat->M_, slumat->colstart[this->cols],
308 slumat->values,slumat->rowindex, slumat->colstart, SLU_NC,
309 static_cast<Dtype_t
>(GetSuperLUType<T>::type), SLU_GE);
312 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:423
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:457
A dense n x m matrix.
Definition: fmatrix.hh:68
Provides access to an iterator over all matrix rows.
Definition: colcompmatrix.hh:22
Provides access to an iterator over an arbitrary subset of matrix rows.
Definition: colcompmatrix.hh:60
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:742
virtual void setMatrix(const Matrix &mat)
Initialize data from given matrix.
Definition: supermatrix.hh:267
virtual ~SuperLUMatrix()
Destructor.
Definition: supermatrix.hh:211
SuperLUMatrix(const Matrix &mat)
Constructor that initializes the data.
Definition: supermatrix.hh:204
BCRSMatrix< FieldMatrix< B, n, m >, TA > Matrix
The type of the matrix to convert.
Definition: supermatrix.hh:194
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:255
virtual void free()
free allocated space.
Definition: supermatrix.hh:277
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:10
Utility class for converting an ISTL Matrix into a column-compressed matrix.
Definition: colcompmatrix.hh:145
Utility class for converting an ISTL Matrix into a SuperLU Matrix.
Definition: supermatrix.hh:173
Traits for type conversions and type information.