HermesCommon  2.0
precond_ml.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_ML
20 #include "precond_ml.h"
21 
22 namespace Hermes
23 {
24  namespace Preconditioners
25  {
26  template<typename Scalar>
27  MlPrecond<Scalar>::MlPrecond(const char *type)
28  {
29  this->prec = NULL;
30  this->owner = true;
31  this->mat = NULL;
32 
33  if(strcmp(type, "sa") == 0) ML_Epetra::SetDefaults("SA", mlist);
34  else if(strcmp(type, "dd") == 0) ML_Epetra::SetDefaults("DD", mlist);
35  }
36 
37  template<typename Scalar>
38  MlPrecond<Scalar>::MlPrecond(ML_Epetra::MultiLevelPreconditioner *mpc)
39  {
40  this->prec = mpc;
41  this->owner = false;
42  this->mat = NULL; // FIXME: get the matrix from mpc
43  }
44 
45  template<typename Scalar>
46  MlPrecond<Scalar>::~MlPrecond()
47  {
48  if(owner) delete prec;
49  }
50 
51  template<typename Scalar>
52  void MlPrecond<Scalar>::set_param(const char *name, const char *value)
53  {
54  mlist.set(name, value);
55  }
56 
57  template<typename Scalar>
58  void MlPrecond<Scalar>::set_param(const char *name, int value)
59  {
60  mlist.set(name, value);
61  }
62 
63  template<typename Scalar>
64  void MlPrecond<Scalar>::set_param(const char *name, double value)
65  {
66  mlist.set(name, value);
67  }
68 
69  template<typename Scalar>
70  void MlPrecond<Scalar>::create(Matrix<Scalar> *m)
71  {
72  EpetraMatrix<Scalar> *mt = static_cast<EpetraMatrix<Scalar> *>(m);
73  assert(mt != NULL);
74  mat = mt;
75  delete prec;
76  prec = new ML_Epetra::MultiLevelPreconditioner(*mat->mat, mlist, false);
77  }
78 
79  template<typename Scalar>
80  void MlPrecond<Scalar>::destroy()
81  {
82  assert(prec != NULL);
83  prec->DestroyPreconditioner();
84  }
85 
86  template<typename Scalar>
87  void MlPrecond<Scalar>::compute()
88  {
89  assert(prec != NULL);
90  prec->ComputePreconditioner();
91  }
92 
93  template<typename Scalar>
94  void MlPrecond<Scalar>::print_unused()
95  {
96  assert(prec != NULL);
97  prec->PrintUnused();
98  }
99 
100  template<typename Scalar>
101  int MlPrecond<Scalar>::ApplyInverse(const Epetra_MultiVector &r, Epetra_MultiVector &z) const
102  {
103  assert(prec != NULL);
104  return prec->ApplyInverse(r, z);
105  }
106 
107  template<typename Scalar>
108  const Epetra_Comm &MlPrecond<Scalar>::Comm() const
109  {
110  return mat->mat->Comm();
111  }
112 
113  template<typename Scalar>
114  const Epetra_Map &MlPrecond<Scalar>::OperatorDomainMap() const
115  {
116  return mat->mat->OperatorDomainMap();
117  }
118 
119  template<typename Scalar>
120  const Epetra_Map &MlPrecond<Scalar>::OperatorRangeMap() const
121  {
122  return mat->mat->OperatorRangeMap();
123  }
124  template class HERMES_API MlPrecond<double>;
125  template class HERMES_API MlPrecond<std::complex<double> >;
126  }
127 }
128 #endif