Dune Core Modules (2.6.0)

ilusubdomainsolver.hh
Go to the documentation of this file.
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_ILUSUBDOMAINSOLVER_HH
4 #define DUNE_ISTL_ILUSUBDOMAINSOLVER_HH
5 
6 #include <map>
8 #include "matrix.hh"
9 #include <cmath>
10 #include <cstdlib>
11 
12 namespace Dune {
13 
32  template<class M, class X, class Y>
34  public:
36  typedef typename std::remove_const<M>::type matrix_type;
38  typedef X domain_type;
40  typedef Y range_type;
41 
48  virtual void apply (X& v, const Y& d) =0;
49 
50  virtual ~ILUSubdomainSolver()
51  {}
52 
53  protected:
59  template<class S>
60  std::size_t copyToLocalMatrix(const M& A, S& rowset);
61 
63  // for ILUN
65  };
66 
73  template<class M, class X, class Y>
75  : public ILUSubdomainSolver<M,X,Y>{
76  public:
78  typedef typename std::remove_const<M>::type matrix_type;
79  typedef typename std::remove_const<M>::type rilu_type;
81  typedef X domain_type;
83  typedef Y range_type;
84 
85 
90  void apply (X& v, const Y& d)
91  {
92  bilu_backsolve(this->ILU,v,d);
93  }
101  template<class S>
102  void setSubMatrix(const M& A, S& rowset);
103 
104  };
105 
106  template<class M, class X, class Y>
107  class ILUNSubdomainSolver
108  : public ILUSubdomainSolver<M,X,Y>{
109  public:
111  typedef typename std::remove_const<M>::type matrix_type;
112  typedef typename std::remove_const<M>::type rilu_type;
114  typedef X domain_type;
116  typedef Y range_type;
117 
122  void apply (X& v, const Y& d)
123  {
124  bilu_backsolve(RILU,v,d);
125  }
126 
134  template<class S>
135  void setSubMatrix(const M& A, S& rowset);
136 
137  private:
141  rilu_type RILU;
142  };
143 
144 
145 
146  template<class M, class X, class Y>
147  template<class S>
148  std::size_t ILUSubdomainSolver<M,X,Y>::copyToLocalMatrix(const M& A, S& rowSet)
149  {
150  // Calculate consecutive indices for local problem
151  // while perserving the ordering
152  typedef typename M::size_type size_type;
153  typedef std::map<typename S::value_type,size_type> IndexMap;
154  typedef typename IndexMap::iterator IMIter;
155  IndexMap indexMap;
156  IMIter guess = indexMap.begin();
157  size_type localIndex=0;
158 
159  typedef typename S::const_iterator SIter;
160  for(SIter rowIdx = rowSet.begin(), rowEnd=rowSet.end();
161  rowIdx!= rowEnd; ++rowIdx, ++localIndex)
162  guess = indexMap.insert(guess,
163  std::make_pair(*rowIdx,localIndex));
164 
165 
166  // Build Matrix for local subproblem
167  ILU.setSize(rowSet.size(),rowSet.size());
168  ILU.setBuildMode(matrix_type::row_wise);
169 
170  // Create sparsity pattern
171  typedef typename matrix_type::CreateIterator CIter;
172  CIter rowCreator = ILU.createbegin();
173  std::size_t offset=0;
174  for(SIter rowIdx = rowSet.begin(), rowEnd=rowSet.end();
175  rowIdx!= rowEnd; ++rowIdx, ++rowCreator) {
176  // See which row entries are in our subset and add them to
177  // the sparsity pattern
178  guess = indexMap.begin();
179 
180  for(typename matrix_type::ConstColIterator col=A[*rowIdx].begin(),
181  endcol=A[*rowIdx].end(); col != endcol; ++col) {
182  // search for the entry in the row set
183  guess = indexMap.find(col.index());
184  if(guess!=indexMap.end()) {
185  // add local index to row
186  rowCreator.insert(guess->second);
187  offset=std::max(offset,(std::size_t)std::abs((int)(guess->second-rowCreator.index())));
188  }
189  }
190 
191  }
192 
193  // Insert the matrix values for the local problem
194  typename matrix_type::iterator iluRow=ILU.begin();
195 
196  for(SIter rowIdx = rowSet.begin(), rowEnd=rowSet.end();
197  rowIdx!= rowEnd; ++rowIdx, ++iluRow) {
198  // See which row entries are in our subset and add them to
199  // the sparsity pattern
200  typename matrix_type::ColIterator localCol=iluRow->begin();
201  for(typename matrix_type::ConstColIterator col=A[*rowIdx].begin(),
202  endcol=A[*rowIdx].end(); col != endcol; ++col) {
203  // search for the entry in the row set
204  guess = indexMap.find(col.index());
205  if(guess!=indexMap.end()) {
206  // set local value
207  (*localCol)=(*col);
208  ++localCol;
209  }
210  }
211  }
212  return offset;
213  }
214 
215 
216  template<class M, class X, class Y>
217  template<class S>
218  void ILU0SubdomainSolver<M,X,Y>::setSubMatrix(const M& A, S& rowSet)
219  {
220  this->copyToLocalMatrix(A,rowSet);
221  bilu0_decomposition(this->ILU);
222  }
223 
224  template<class M, class X, class Y>
225  template<class S>
226  void ILUNSubdomainSolver<M,X,Y>::setSubMatrix(const M& A, S& rowSet)
227  {
228  std::size_t offset=copyToLocalMatrix(A,rowSet);
229  RILU.setSize(rowSet.size(),rowSet.size(), (1+2*offset)*rowSet.size());
230  RILU.setBuildMode(matrix_type::row_wise);
231  bilu_decomposition(this->ILU, (offset+1)/2, RILU);
232  }
233 
235 } // end name space DUNE
236 
237 
238 #endif
Exact subdomain solver using ILU(p) with appropriate p.
Definition: ilusubdomainsolver.hh:75
X domain_type
The domain type of the preconditioner.
Definition: ilusubdomainsolver.hh:81
Y range_type
The range type of the preconditioner.
Definition: ilusubdomainsolver.hh:83
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: ilusubdomainsolver.hh:78
void apply(X &v, const Y &d)
Apply the subdomain solver.
Definition: ilusubdomainsolver.hh:90
base class encapsulating common algorithms of ILU0SubdomainSolver and ILUNSubdomainSolver.
Definition: ilusubdomainsolver.hh:33
matrix_type ILU
The ILU0 decomposition of the matrix, or the local matrix.
Definition: ilusubdomainsolver.hh:64
X domain_type
The domain type of the preconditioner.
Definition: ilusubdomainsolver.hh:38
Y range_type
The range type of the preconditioner.
Definition: ilusubdomainsolver.hh:40
std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: ilusubdomainsolver.hh:36
virtual void apply(X &v, const Y &d)=0
Apply the subdomain solver.
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:58
void bilu_backsolve(const M &A, X &v, const Y &d)
LU backsolve with stored inverse.
Definition: ilu.hh:97
void bilu_decomposition(const M &A, int n, M &ILU)
Definition: ilu.hh:157
void bilu0_decomposition(M &A)
compute ILU decomposition of A. A is overwritten by its decomposition
Definition: ilu.hh:36
std::size_t copyToLocalMatrix(const M &A, S &rowset)
Copy the local part of the global matrix to ILU.
Definition: ilusubdomainsolver.hh:148
void setSubMatrix(const M &A, S &rowset)
Set the data of the local problem.
Definition: ilusubdomainsolver.hh:218
A dynamic dense block matrix class.
Dune namespace.
Definition: alignedallocator.hh:10
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 2, 22:35, 2024)