HermesCommon  2.0
precond_ifpack.cpp
Go to the documentation of this file.
1 // This file is part of HermesCommon
2 //
3 // Hermes2D is free software: you can redistribute it and/or modify
4 // it under the terms of the GNU General Public License as published by
5 // the Free Software Foundation, either version 2 of the License, or
6 // (at your option) any later version.
7 //
8 // Hermes2D is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11 // GNU General Public License for more details.
12 //
13 // You should have received a copy of the GNU General Public License
14 // along with Hermes2D. If not, see <http://www.gnu.org/licenses/>.
18 #include "config.h"
19 #ifdef HAVE_IFPACK
20 #include "precond_ifpack.h"
21 #include "Ifpack_PointRelaxation.h"
22 #include "Ifpack_BlockRelaxation.h"
23 #include "Ifpack_DenseContainer.h"
24 #include "Ifpack_AdditiveSchwarz.h"
25 #include "Ifpack_ILU.h"
26 #include "Ifpack_ILUT.h"
27 #include "Ifpack_IC.h"
28 #include "Ifpack_ICT.h"
29 #include "Ifpack_Graph_Epetra_CrsGraph.h"
30 
31 namespace Hermes
32 {
33  namespace Preconditioners
34 {
35  template<typename Scalar>
36  IfpackPrecond<Scalar>::IfpackPrecond(const char *cls, const char *type)
37  {
38  this->prec = NULL;
39  this->owner = true;
40  this->mat = NULL;
41 
42  this->cls = cls;
43  this->type = type;
44  }
45 
46  template<typename Scalar>
47  IfpackPrecond<Scalar>::IfpackPrecond(const char *cls, const char *type, int overlap)
48  {
49  this->prec = NULL;
50  this->owner = true;
51  this->mat = NULL;
52 
53  this->cls = cls;
54  this->type = type;
55  this->overlap = overlap;
56  }
57 
58  template<typename Scalar>
59  IfpackPrecond<Scalar>::IfpackPrecond(Ifpack_Preconditioner *ipc)
60  {
61  this->prec = ipc;
62  this->owner = false;
63  this->mat = NULL; // FIXME: take the matrix from ipc
64  }
65 
66  template<typename Scalar>
67  IfpackPrecond<Scalar>::~IfpackPrecond()
68  {
69  if(owner) delete prec;
70  }
71 
72  template<typename Scalar>
73  void IfpackPrecond<Scalar>::set_param(const char *name, const char *value)
74  {
75  ilist.set(name, value);
76  }
77 
78  template<typename Scalar>
79  void IfpackPrecond<Scalar>::set_param(const char *name, int value)
80  {
81  ilist.set(name, value);
82  }
83 
84  template<typename Scalar>
85  void IfpackPrecond<Scalar>::set_param(const char *name, double value)
86  {
87  ilist.set(name, value);
88  }
89 
90  template<typename Scalar>
91  void IfpackPrecond<Scalar>::create(Matrix<Scalar> *m)
92  {
93  EpetraMatrix<Scalar> *mt = static_cast<EpetraMatrix<Scalar> *>(m);
94  assert(mt != NULL);
95  mat = mt;
96  if(strcmp(cls, "point-relax") == 0)
97  {
98  create_point_relax(mat, type);
99  apply_params();
100  initialize();
101  }
102  else if(strcmp(cls, "block-relax") == 0)
103  {
104  create_block_relax(mat, type);
105  apply_params();
106  }
107  else if(strcmp(cls, "add-schwartz") == 0)
108  {
109  create_add_schwartz(mat, type, overlap);
110  apply_params();
111  initialize();
112  }
113  }
114 
115  template<typename Scalar>
116  void IfpackPrecond<Scalar>::create_point_relax(EpetraMatrix<Scalar> *a, const char *name)
117  {
118  prec = new Ifpack_PointRelaxation(a->mat);
119  ilist.set("relaxation: type", name);
120  }
121 
122  template<typename Scalar>
123  void IfpackPrecond<Scalar>::create_block_relax(EpetraMatrix<Scalar> *a, const char *name)
124  {
125  Teuchos::RCP<const Epetra_CrsGraph> rgraph = Teuchos::rcp(&(a->mat->Graph()));
126  Ifpack_Graph *graph = new Ifpack_Graph_Epetra_CrsGraph(rgraph);
127 
128  Ifpack_Partitioner *partitioner = new Ifpack_GreedyPartitioner(graph);
129 
130  Teuchos::ParameterList list;
131  list.set("partitioner: local parts", 1000); //\todo parametrize me
132  partitioner->SetParameters(list);
133  partitioner->Compute();
134 
135  prec = new Ifpack_BlockRelaxation<Ifpack_DenseContainer>(a->mat);
136  ilist.set("relaxation: type", name);
137 
138  rgraph.release();
139  }
140 
141  template<typename Scalar>
142  void IfpackPrecond<Scalar>::create_add_schwartz(EpetraMatrix<Scalar> *a, const char *name, int overlap)
143  {
144  if(strcasecmp(name, "ilu") == 0)
145  {
146  prec = new Ifpack_AdditiveSchwarz<Ifpack_ILU>(a->mat, overlap);
147  }
148  else if(strcasecmp(name, "ilut") == 0)
149  {
150  prec = new Ifpack_AdditiveSchwarz<Ifpack_ILUT>(a->mat, overlap);
151  }
152  else if(strcasecmp(name, "ic") == 0)
153  {
154  prec = new Ifpack_AdditiveSchwarz<Ifpack_IC>(a->mat, overlap);
155  }
156  else if(strcasecmp(name, "ict") == 0)
157  {
158  prec = new Ifpack_AdditiveSchwarz<Ifpack_ICT>(a->mat, overlap);
159  }
160  else
161  prec = NULL;
162  }
163 
164  template<typename Scalar>
165  int IfpackPrecond<Scalar>::initialize()
166  {
167  assert(prec != NULL);
168  return prec->Initialize();
169  }
170 
171  template<typename Scalar>
172  void IfpackPrecond<Scalar>::compute()
173  {
174  assert(prec != NULL);
175  prec->Compute();
176  }
177 
178  template<typename Scalar>
179  void IfpackPrecond<Scalar>::apply_params()
180  {
181  prec->SetParameters(ilist);
182  }
183 
184  template<typename Scalar>
185  int IfpackPrecond<Scalar>::ApplyInverse(const Epetra_MultiVector &r, Epetra_MultiVector &z) const
186  {
187  assert(prec != NULL);
188  return prec->ApplyInverse(r, z);
189  }
190 
191  template<typename Scalar>
192  const Epetra_Comm &IfpackPrecond<Scalar>::Comm() const
193  {
194  return mat->mat->Comm();
195  }
196 
197  template<typename Scalar>
198  const Epetra_Map &IfpackPrecond<Scalar>::OperatorDomainMap() const
199  {
200  return mat->mat->OperatorDomainMap();
201  }
202 
203  template<typename Scalar>
204  const Epetra_Map &IfpackPrecond<Scalar>::OperatorRangeMap() const
205  {
206  return mat->mat->OperatorRangeMap();
207  }
208  template class HERMES_API IfpackPrecond<double>;
209  template class HERMES_API IfpackPrecond<std::complex<double> >;
210  }
211 }
212 #endif