dune-istl  2.2.1
ilu.hh
Go to the documentation of this file.
1 #ifndef DUNE_ILU_HH
2 #define DUNE_ILU_HH
3 
4 #include<cmath>
5 #include<complex>
6 #include<iostream>
7 #include<iomanip>
8 #include<string>
9 #include<set>
10 #include<map>
11 
12 #include "istlexception.hh"
13 #include "io.hh"
14 
19 namespace Dune {
20 
25  class MatrixBlockError : public virtual Dune::FMatrixError {
26  public:
27  int r, c;
28  };
29 
30 
32  template<class M>
34  {
35  // iterator types
36  typedef typename M::RowIterator rowiterator;
37  typedef typename M::ColIterator coliterator;
38  typedef typename M::block_type block;
39 
40  // implement left looking variant with stored inverse
41  rowiterator endi=A.end();
42  for (rowiterator i=A.begin(); i!=endi; ++i)
43  {
44  // coliterator is diagonal after the following loop
45  coliterator endij=(*i).end(); // end of row i
46  coliterator ij;
47 
48  // eliminate entries left of diagonal; store L factor
49  for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
50  {
51  // find A_jj which eliminates A_ij
52  coliterator jj = A[ij.index()].find(ij.index());
53 
54  // compute L_ij = A_jj^-1 * A_ij
55  (*ij).rightmultiply(*jj);
56 
57  // modify row
58  coliterator endjk=A[ij.index()].end(); // end of row j
59  coliterator jk=jj; ++jk;
60  coliterator ik=ij; ++ik;
61  while (ik!=endij && jk!=endjk)
62  if (ik.index()==jk.index())
63  {
64  block B(*jk);
65  B.leftmultiply(*ij);
66  *ik -= B;
67  ++ik; ++jk;
68  }
69  else
70  {
71  if (ik.index()<jk.index())
72  ++ik;
73  else
74  ++jk;
75  }
76  }
77 
78  // invert pivot and store it in A
79  if (ij.index()!=i.index())
80  DUNE_THROW(ISTLError,"diagonal entry missing");
81  try {
82  (*ij).invert(); // compute inverse of diagonal block
83  }
84  catch (Dune::FMatrixError & e) {
85  DUNE_THROW(MatrixBlockError, "ILU failed to invert matrix block A["
86  << i.index() << "][" << ij.index() << "]" << e.what();
87  th__ex.r=i.index(); th__ex.c=ij.index(););
88  }
89  }
90  }
91 
93  template<class M, class X, class Y>
94  void bilu_backsolve (const M& A, X& v, const Y& d)
95  {
96  // iterator types
97  typedef typename M::ConstRowIterator rowiterator;
98  typedef typename M::ConstColIterator coliterator;
99  typedef typename Y::block_type dblock;
100  typedef typename X::block_type vblock;
101 
102  // lower triangular solve
103  rowiterator endi=A.end();
104  for (rowiterator i=A.begin(); i!=endi; ++i)
105  {
106  dblock rhs(d[i.index()]);
107  for (coliterator j=(*i).begin(); j.index()<i.index(); ++j)
108  (*j).mmv(v[j.index()],rhs);
109  v[i.index()] = rhs; // Lii = I
110  }
111 
112  // upper triangular solve
113  rowiterator rendi=A.beforeBegin();
114  for (rowiterator i=A.beforeEnd(); i!=rendi; --i)
115  {
116  vblock rhs(v[i.index()]);
117  coliterator j;
118  for (j=(*i).beforeEnd(); j.index()>i.index(); --j)
119  (*j).mmv(v[j.index()],rhs);
120  v[i.index()] = 0;
121  (*j).umv(rhs,v[i.index()]); // diagonal stores inverse!
122  }
123  }
124 
125 
126 
127  // recursive function template to access first entry of a matrix
128  template<class M>
129  typename M::field_type& firstmatrixelement (M& A)
130  {
131  return firstmatrixelement(*(A.begin()->begin()));
132  }
133 
134  template<class K, int n, int m>
135  K& firstmatrixelement (FieldMatrix<K,n,m>& A)
136  {
137  return A[0][0];
138  }
139 
140  template<class K>
141  K& firstmatrixelement (FieldMatrix<K,1,1>& A)
142  {
143  return A[0][0];
144  }
145 
146 
153  template<class M>
154  void bilu_decomposition (const M& A, int n, M& ILU)
155  {
156  // iterator types
157  typedef typename M::RowIterator rowiterator;
158  typedef typename M::ColIterator coliterator;
159  typedef typename M::ConstRowIterator crowiterator;
160  typedef typename M::ConstColIterator ccoliterator;
161  typedef typename M::CreateIterator createiterator;
162  typedef typename M::block_type block;
163  typedef typename M::field_type K;
164  typedef std::map<size_t, int> map;
165  typedef typename map::iterator mapiterator;
166 
167  // symbolic factorization phase, store generation number in first matrix element
168  crowiterator endi=A.end();
169  createiterator ci=ILU.createbegin();
170  for (crowiterator i=A.begin(); i!=endi; ++i)
171  {
172  // std::cout << "in row " << i.index() << std::endl;
173  map rowpattern; // maps column index to generation
174 
175  // initialize pattern with row of A
176  for (ccoliterator j=(*i).begin(); j!=(*i).end(); ++j)
177  rowpattern[j.index()] = 0;
178 
179  // eliminate entries in row which are to the left of the diagonal
180  for (mapiterator ik=rowpattern.begin(); (*ik).first<i.index(); ++ik)
181  {
182  if ((*ik).second<n)
183  {
184 // std::cout << " eliminating " << i.index() << "," << (*ik).first
185 // << " level " << (*ik).second << std::endl;
186 
187  coliterator endk = ILU[(*ik).first].end(); // end of row k
188  coliterator kj = ILU[(*ik).first].find((*ik).first); // diagonal in k
189  for (++kj; kj!=endk; ++kj) // row k eliminates in row i
190  {
191  int generation = (int) firstmatrixelement(*kj);
192  if (generation<n)
193  {
194  mapiterator ij = rowpattern.find(kj.index());
195  if (ij==rowpattern.end())
196  {
197  //std::cout << " new entry " << i.index() << "," << kj.index()
198 // << " level " << (*ik).second+1 << std::endl;
199 
200  rowpattern[kj.index()] = generation+1;
201  }
202  }
203  }
204  }
205  }
206 
207  // create row
208  for (mapiterator ik=rowpattern.begin(); ik!=rowpattern.end(); ++ik)
209  ci.insert((*ik).first);
210  ++ci; // now row i exist
211 
212  // write generation index into entries
213  coliterator endILUij = ILU[i.index()].end();;
214  for (coliterator ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
215  firstmatrixelement(*ILUij) = (K) rowpattern[ILUij.index()];
216  }
217 
218  // printmatrix(std::cout,ILU,"ilu pattern","row",10,2);
219 
220  // copy entries of A
221  for (crowiterator i=A.begin(); i!=endi; ++i)
222  {
223  coliterator ILUij;
224  coliterator endILUij = ILU[i.index()].end();;
225  for (ILUij=ILU[i.index()].begin(); ILUij!=endILUij; ++ILUij)
226  (*ILUij) = 0; // clear row
227  ccoliterator Aij = (*i).begin();
228  ccoliterator endAij = (*i).end();
229  ILUij = ILU[i.index()].begin();
230  while (Aij!=endAij && ILUij!=endILUij)
231  {
232  if (Aij.index()==ILUij.index())
233  {
234  *ILUij = *Aij;
235  ++Aij; ++ILUij;
236  }
237  else
238  {
239  if (Aij.index()<ILUij.index())
240  ++Aij;
241  else
242  ++ILUij;
243  }
244  }
245  }
246 
247  // call decomposition on pattern
248  bilu0_decomposition(ILU);
249  }
250 
251 
254 } // end namespace
255 
256 #endif