Dune Core Modules (2.4.1)

supermatrix.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ISTL_SUPERMATRIX_HH
4#define DUNE_ISTL_SUPERMATRIX_HH
5
6#if HAVE_SUPERLU
7#ifdef SUPERLU_POST_2005_VERSION
8
9#ifndef SUPERLU_NTYPE
10#define SUPERLU_NTYPE 1
11#endif
12
13#if SUPERLU_NTYPE==0
14#include "slu_sdefs.h"
15#endif
16
17#if SUPERLU_NTYPE==1
18#include "slu_ddefs.h"
19#endif
20
21#if SUPERLU_NTYPE==2
22#include "slu_cdefs.h"
23#endif
24
25#if SUPERLU_NTYPE>=3
26#include "slu_zdefs.h"
27#endif
28
29#else
30
31#if SUPERLU_NTYPE==0
32#include "ssp_defs.h"
33#endif
34
35#if SUPERLU_NTYPE==1
36#include "dsp_defs.h"
37#endif
38
39#if SUPERLU_NTYPE==2
40#include "csp_defs.h"
41#endif
42
43#if SUPERLU_NTYPE>=3
44#include "zsp_defs.h"
45#endif
46
47#endif
48#include "bcrsmatrix.hh"
49#include "bvector.hh"
53#include <limits>
54
55#include"colcompmatrix.hh"
56
57namespace Dune
58{
59
60 template<class T>
61 struct SuperMatrixCreateSparseChooser
62 {};
63
64 template<class T>
65 struct SuperMatrixPrinter
66 {};
67
68#if SUPERLU_NTYPE==0
69 template<>
70 struct SuperMatrixCreateSparseChooser<float>
71 {
72 static void create(SuperMatrix *mat, int n, int m, int offset,
73 float *values, int *rowindex, int* colindex,
74 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
75 {
76 sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
77 stype, dtype, mtype);
78 }
79 };
80
81 template<>
82 struct SuperMatrixPrinter<float>
83 {
84 static void print(char* name, SuperMatrix* mat)
85 {
86 sPrint_CompCol_Matrix(name, mat);
87 }
88 };
89#endif
90
91#if SUPERLU_NTYPE==1
92 template<>
93 struct SuperMatrixCreateSparseChooser<double>
94 {
95 static void create(SuperMatrix *mat, int n, int m, int offset,
96 double *values, int *rowindex, int* colindex,
97 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
98 {
99 dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
100 stype, dtype, mtype);
101 }
102 };
103
104 template<>
105 struct SuperMatrixPrinter<double>
106 {
107 static void print(char* name, SuperMatrix* mat)
108 {
109 dPrint_CompCol_Matrix(name, mat);
110 }
111 };
112#endif
113
114#if SUPERLU_NTYPE==2
115 template<>
116 struct SuperMatrixCreateSparseChooser<std::complex<float> >
117 {
118 static void create(SuperMatrix *mat, int n, int m, int offset,
119 std::complex<float> *values, int *rowindex, int* colindex,
120 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
121 {
122 cCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast< ::complex*>(values),
123 rowindex, colindex, stype, dtype, mtype);
124 }
125 };
126
127 template<>
128 struct SuperMatrixPrinter<std::complex<float> >
129 {
130 static void print(char* name, SuperMatrix* mat)
131 {
132 cPrint_CompCol_Matrix(name, mat);
133 }
134 };
135#endif
136
137#if SUPERLU_NTYPE>=3
138 template<>
139 struct SuperMatrixCreateSparseChooser<std::complex<double> >
140 {
141 static void create(SuperMatrix *mat, int n, int m, int offset,
142 std::complex<double> *values, int *rowindex, int* colindex,
143 Stype_t stype, Dtype_t dtype, Mtype_t mtype)
144 {
145 zCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast<doublecomplex*>(values),
146 rowindex, colindex, stype, dtype, mtype);
147 }
148 };
149
150 template<>
151 struct SuperMatrixPrinter<std::complex<double> >
152 {
153 static void print(char* name, SuperMatrix* mat)
154 {
155 zPrint_CompCol_Matrix(name, mat);
156 }
157 };
158#endif
159
160 template<class T>
161 struct BaseGetSuperLUType
162 {
163 static const Dtype_t type;
164 };
165
166 template<class T>
167 struct GetSuperLUType
168 {};
169
170 template<class T>
171 const Dtype_t BaseGetSuperLUType<T>::type =
172 Dune::is_same<T,float>::value ? SLU_S :
173 ( Dune::is_same<T,std::complex<double> >::value ? SLU_Z :
174 ( Dune::is_same<T,std::complex<float> >::value ? SLU_C : SLU_D ));
175
176 template<>
177 struct GetSuperLUType<double>
178 : public BaseGetSuperLUType<double>
179 {
180 typedef double float_type;
181 };
182
183 template<>
184 struct GetSuperLUType<float>
185 : public BaseGetSuperLUType<float>
186 {
187 typedef float float_type;
188 };
189
190 template<>
191 struct GetSuperLUType<std::complex<double> >
192 : public BaseGetSuperLUType<std::complex<double> >
193 {
194 typedef double float_type;
195 };
196
197 template<>
198 struct GetSuperLUType<std::complex<float> >
199 : public BaseGetSuperLUType<std::complex<float> >
200 {
201 typedef float float_type;
202
203 };
204
209 template<class M>
211 {};
212
213 template<class M>
214 struct SuperMatrixInitializer
215 {};
216
217 template<class T>
218 class SuperLU;
219
223 template<class B, class TA, int n, int m>
225 : public ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
226 {
227 template<class M, class X, class TM, class TD, class T1>
228 friend class SeqOverlappingSchwarz;
229 friend struct SuperMatrixInitializer<BCRSMatrix<FieldMatrix<B,n,m>,TA> >;
230 public:
233
234 friend struct SeqOverlappingSchwarzAssembler<SuperLU<Matrix> >;
235
236 typedef typename Matrix::size_type size_type;
237
242 explicit SuperLUMatrix(const Matrix& mat) : ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >(mat)
243 {}
244
246 {}
247
250 {
251 if (this->N_+this->M_*this->Nnz_ != 0)
252 free();
253 }
254
256 operator SuperMatrix&()
257 {
258 return A;
259 }
260
262 operator const SuperMatrix&() const
263 {
264 return A;
265 }
266
268 {
269 this->ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(mat);
270 SuperMatrixCreateSparseChooser<B>
271 ::create(&A, this->N_, this->M_, this->colstart[this->N_],
272 this->values,this->rowindex, this->colstart, SLU_NC,
273 static_cast<Dtype_t>(GetSuperLUType<B>::type), SLU_GE);
274 return *this;
275 }
276
277 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& operator=(const SuperLUMatrix <BCRSMatrix<FieldMatrix<B,n,m>,TA> >& mat)
278 {
279 this->ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(mat);
280 SuperMatrixCreateSparseChooser<B>
281 ::create(&A, this->N_, this->M_, this->colstart[this->N_],
282 this->values,this->rowindex, this->colstart, SLU_NC,
283 static_cast<Dtype_t>(GetSuperLUType<B>::type), SLU_GE);
284 return *this;
285 }
286
293 virtual void setMatrix(const Matrix& mat, const std::set<std::size_t>& mrs)
294 {
295 if(this->N_+this->M_+this->Nnz_!=0)
296 free();
297 this->N_=mrs.size()*n;
298 this->M_=mrs.size()*m;
299 SuperMatrixInitializer<Matrix> initializer(*this);
300
301 copyToColCompMatrix(initializer, MatrixRowSubset<Matrix,std::set<std::size_t> >(mat,mrs));
302 }
303
305 virtual void setMatrix(const Matrix& mat)
306 {
307 this->N_=n*mat.N();
308 this->M_=m*mat.M();
309 SuperMatrixInitializer<Matrix> initializer(*this);
310
311 copyToColCompMatrix(initializer, MatrixRowSet<Matrix>(mat));
312 }
313
315 virtual void free()
316 {
318 SUPERLU_FREE(A.Store);
319 }
320 private:
321 SuperMatrix A;
322 };
323
324 template<class T, class A, int n, int m>
325 class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
326 : public ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
327 {
328 template<class I, class S, class D>
329 friend class OverlappingSchwarzInitializer;
330 public:
331 typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
332 typedef Dune::SuperLUMatrix<Matrix> SuperLUMatrix;
333
334 SuperMatrixInitializer(SuperLUMatrix& lum) : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >(lum)
335 ,slumat(&lum)
336 {}
337
338 SuperMatrixInitializer() : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >()
339 {}
340
341 virtual void createMatrix() const
342 {
343 ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::createMatrix();
344 SuperMatrixCreateSparseChooser<T>
345 ::create(&slumat->A, slumat->N_, slumat->M_, slumat->colstart[this->cols],
346 slumat->values,slumat->rowindex, slumat->colstart, SLU_NC,
347 static_cast<Dtype_t>(GetSuperLUType<T>::type), SLU_GE);
348 }
349 private:
350 SuperLUMatrix* slumat;
351 };
352}
353#endif
354#endif
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:413
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:447
A dense n x m matrix.
Definition: fmatrix.hh:67
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:753
virtual void setMatrix(const Matrix &mat)
Initialize data from given matrix.
Definition: supermatrix.hh:305
virtual ~SuperLUMatrix()
Destructor.
Definition: supermatrix.hh:249
SuperLUMatrix(const Matrix &mat)
Constructor that initializes the data.
Definition: supermatrix.hh:242
BCRSMatrix< FieldMatrix< B, n, m >, TA > Matrix
The type of the matrix to convert.
Definition: supermatrix.hh:232
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:293
virtual void free()
free allocated space.
Definition: supermatrix.hh:315
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: alignment.hh:10
STL namespace.
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:211
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)