HermesCommon  3.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 = nullptr;
30  this->owner = true;
31  this->mat = nullptr;
32 
33  if (strcmp(type, "sa") == 0) ML_Epetra::SetDefaults("SA", mlist);
34  else if (strcmp(type, "nssa") == 0) ML_Epetra::SetDefaults("NSSA", mlist);
35  else if (strcmp(type, "dd") == 0) ML_Epetra::SetDefaults("DD", mlist);
36  }
37 
38  template<typename Scalar>
39  MlPrecond<Scalar>::MlPrecond(ML_Epetra::MultiLevelPreconditioner *mpc)
40  {
41  this->prec = mpc;
42  this->owner = false;
43  // FIXME: get the matrix from mpc
44  this->mat = nullptr;
45  }
46 
47  template<typename Scalar>
49  {
50  if (owner) delete prec;
51  }
52 
53  template<typename Scalar>
54  void MlPrecond<Scalar>::set_param(const char *name, const char *value)
55  {
56  mlist.set(name, value);
57  }
58 
59  template<typename Scalar>
60  void MlPrecond<Scalar>::set_param(const char *name, int value)
61  {
62  mlist.set(name, value);
63  }
64 
65  template<typename Scalar>
66  void MlPrecond<Scalar>::set_param(const char *name, bool value)
67  {
68  mlist.set(name, value);
69  }
70 
71  template<typename Scalar>
72  void MlPrecond<Scalar>::set_param(const char *name, double value)
73  {
74  mlist.set(name, value);
75  }
76 
77  template<typename Scalar>
79  {
80  EpetraMatrix<Scalar> *mt = static_cast<EpetraMatrix<Scalar> *>(m);
81  assert(mt != nullptr);
82  mat = mt;
83  if (prec) delete prec;
84  prec = new ML_Epetra::MultiLevelPreconditioner(*mat->mat, mlist, false);
85  }
86 
87  template<typename Scalar>
89  {
90  assert(prec != nullptr);
91  prec->DestroyPreconditioner();
92  }
93 
94  template<typename Scalar>
96  {
97  assert(prec != nullptr);
98  prec->ComputePreconditioner();
99  }
100 
101  template<typename Scalar>
103  {
104  assert(prec != nullptr);
105  assert(prec->IsPreconditionerComputed());
106  prec->ReComputePreconditioner();
107  }
108 
109  template<typename Scalar>
111  {
112  assert(prec != nullptr);
113  prec->PrintUnused();
114  }
115 
116  template<typename Scalar>
117  int MlPrecond<Scalar>::ApplyInverse(const Epetra_MultiVector &r, Epetra_MultiVector &z) const
118  {
119  assert(prec != nullptr);
120  return prec->ApplyInverse(r, z);
121  }
122 
123  template<typename Scalar>
124  const Epetra_Comm &MlPrecond<Scalar>::Comm() const
125  {
126  return mat->mat->Comm();
127  }
128 
129  template<typename Scalar>
130  const Epetra_Map &MlPrecond<Scalar>::OperatorDomainMap() const
131  {
132  return mat->mat->OperatorDomainMap();
133  }
134 
135  template<typename Scalar>
136  const Epetra_Map &MlPrecond<Scalar>::OperatorRangeMap() const
137  {
138  return mat->mat->OperatorRangeMap();
139  }
140  template class HERMES_API MlPrecond < double > ;
141  template class HERMES_API MlPrecond < std::complex<double> > ;
142  }
143 }
144 #endif
General (abstract) matrix representation in Hermes.
Definition: matrix.h:36
General namespace for the Hermes library.
virtual void destroy()
Destroy the preconditioner object.
Definition: precond_ml.cpp:88
virtual void create(Matrix< Scalar > *mat)
Definition: precond_ml.cpp:78
virtual void compute()
Compute the preconditioner.
Definition: precond_ml.cpp:95
Preconditioners built on ML.
Definition: epetra.h:51
virtual void recompute()
Cheaply recompute the preconditioner if matrix values have changed but not their non-zero structure...
Definition: precond_ml.cpp:102
ML (Trilinos package) preconditioners interface.