Dune Core Modules (2.6.0)

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 
8 #include "bcrsmatrix.hh"
9 #include "bvector.hh"
10 #include <dune/common/fmatrix.hh>
11 #include <dune/common/fvector.hh>
13 #include <limits>
14 
15 #include"colcompmatrix.hh"
16 
17 #include "superlufunctions.hh"
18 
19 namespace Dune
20 {
21 
22  template<class T>
23  struct SuperMatrixCreateSparseChooser
24  {};
25 
26  template<class T>
27  struct SuperMatrixPrinter
28  {};
29 
30 #if HAVE_SLU_SDEFS_H
31  template<>
32  struct SuperMatrixCreateSparseChooser<float>
33  {
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)
37  {
38  sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
39  stype, dtype, mtype);
40  }
41  };
42 
43  template<>
44  struct SuperMatrixPrinter<float>
45  {
46  static void print(char* name, SuperMatrix* mat)
47  {
48  sPrint_CompCol_Matrix(name, mat);
49  }
50  };
51 #endif
52 
53 #if HAVE_SLU_DDEFS_H
54  template<>
55  struct SuperMatrixCreateSparseChooser<double>
56  {
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)
60  {
61  dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
62  stype, dtype, mtype);
63  }
64  };
65 
66  template<>
67  struct SuperMatrixPrinter<double>
68  {
69  static void print(char* name, SuperMatrix* mat)
70  {
71  dPrint_CompCol_Matrix(name, mat);
72  }
73  };
74 #endif
75 
76 #if HAVE_SLU_CDEFS_H
77  template<>
78  struct SuperMatrixCreateSparseChooser<std::complex<float> >
79  {
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)
83  {
84  cCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast< ::complex*>(values),
85  rowindex, colindex, stype, dtype, mtype);
86  }
87  };
88 
89  template<>
90  struct SuperMatrixPrinter<std::complex<float> >
91  {
92  static void print(char* name, SuperMatrix* mat)
93  {
94  cPrint_CompCol_Matrix(name, mat);
95  }
96  };
97 #endif
98 
99 #if HAVE_SLU_ZDEFS_H
100  template<>
101  struct SuperMatrixCreateSparseChooser<std::complex<double> >
102  {
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)
106  {
107  zCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast<doublecomplex*>(values),
108  rowindex, colindex, stype, dtype, mtype);
109  }
110  };
111 
112  template<>
113  struct SuperMatrixPrinter<std::complex<double> >
114  {
115  static void print(char* name, SuperMatrix* mat)
116  {
117  zPrint_CompCol_Matrix(name, mat);
118  }
119  };
120 #endif
121 
122  template<class T>
123  struct BaseGetSuperLUType
124  {
125  static const Dtype_t type;
126  };
127 
128  template<class T>
129  struct GetSuperLUType
130  {};
131 
132  template<class T>
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 ));
137 
138  template<>
139  struct GetSuperLUType<double>
140  : public BaseGetSuperLUType<double>
141  {
142  typedef double float_type;
143  };
144 
145  template<>
146  struct GetSuperLUType<float>
147  : public BaseGetSuperLUType<float>
148  {
149  typedef float float_type;
150  };
151 
152  template<>
153  struct GetSuperLUType<std::complex<double> >
154  : public BaseGetSuperLUType<std::complex<double> >
155  {
156  typedef double float_type;
157  };
158 
159  template<>
160  struct GetSuperLUType<std::complex<float> >
161  : public BaseGetSuperLUType<std::complex<float> >
162  {
163  typedef float float_type;
164 
165  };
166 
171  template<class M>
173  {};
174 
175  template<class M>
176  struct SuperMatrixInitializer
177  {};
178 
179  template<class T>
180  class SuperLU;
181 
185  template<class B, class TA, int n, int m>
187  : public ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
188  {
189  template<class M, class X, class TM, class TD, class T1>
190  friend class SeqOverlappingSchwarz;
191  friend struct SuperMatrixInitializer<BCRSMatrix<FieldMatrix<B,n,m>,TA> >;
192  public:
195 
196  friend struct SeqOverlappingSchwarzAssemblerHelper<SuperLU<Matrix>, true>;
197 
198  typedef typename Matrix::size_type size_type;
199 
204  explicit SuperLUMatrix(const Matrix& mat) : ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >(mat)
205  {}
206 
208  {}
209 
211  virtual ~SuperLUMatrix()
212  {
213  if (this->N_+this->M_*this->Nnz_ != 0)
214  free();
215  }
216 
218  operator SuperMatrix&()
219  {
220  return A;
221  }
222 
224  operator const SuperMatrix&() const
225  {
226  return A;
227  }
228 
230  {
231  this->ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(mat);
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);
236  return *this;
237  }
238 
239  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& operator=(const SuperLUMatrix <BCRSMatrix<FieldMatrix<B,n,m>,TA> >& mat)
240  {
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);
246  return *this;
247  }
248 
255  virtual void setMatrix(const Matrix& mat, const std::set<std::size_t>& mrs)
256  {
257  if(this->N_+this->M_+this->Nnz_!=0)
258  free();
259  this->N_=mrs.size()*n;
260  this->M_=mrs.size()*m;
261  SuperMatrixInitializer<Matrix> initializer(*this);
262 
263  copyToColCompMatrix(initializer, MatrixRowSubset<Matrix,std::set<std::size_t> >(mat,mrs));
264  }
265 
267  virtual void setMatrix(const Matrix& mat)
268  {
269  this->N_=n*mat.N();
270  this->M_=m*mat.M();
271  SuperMatrixInitializer<Matrix> initializer(*this);
272 
273  copyToColCompMatrix(initializer, MatrixRowSet<Matrix>(mat));
274  }
275 
277  virtual void free()
278  {
280  SUPERLU_FREE(A.Store);
281  }
282  private:
283  SuperMatrix A;
284  };
285 
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> >
289  {
290  template<class I, class S, class D>
291  friend class OverlappingSchwarzInitializer;
292  public:
293  typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
294  typedef Dune::SuperLUMatrix<Matrix> SuperLUMatrix;
295 
296  SuperMatrixInitializer(SuperLUMatrix& lum) : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >(lum)
297  ,slumat(&lum)
298  {}
299 
300  SuperMatrixInitializer() : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >()
301  {}
302 
303  virtual void createMatrix() const
304  {
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);
310  }
311  private:
312  SuperLUMatrix* slumat;
313  };
314 }
315 #endif // HAVE_SUPERLU
316 #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: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.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 3, 22:32, 2024)