AMDiS  0.1
The Adaptive Multi-Dimensional Simulation Toolbox
PetscRunner.hpp
1 #pragma once
2 
3 #include <petscpc.h>
4 #include <petscksp.h>
5 
6 #include <dune/common/unused.hh>
7 
8 #include <amdis/Initfile.hpp>
9 #include <amdis/linearalgebra/RunnerInterface.hpp>
10 #include <amdis/linearalgebra/SolverInfo.hpp>
11 #include <amdis/linearalgebra/petsc/Interface.hpp>
12 
13 #ifndef KSPDGGMRES
14 #define KSPDGGMRES "dggmres"
15 #endif
16 
17 
18 namespace AMDiS
19 {
21 
45  template <class Matrix, class Vector>
47  : public RunnerInterface<Matrix,Vector>
48  {
49  using M = typename Matrix::BaseMatrix;
50  using X = typename Vector::BaseVector;
51  using Y = typename Vector::BaseVector;
52 
53  public:
55 
63  explicit PetscRunner(std::string const& prefix)
64  : prefix_(prefix)
65  {
66  Parameters::get(prefix + "->info", info_);
67  }
68 
69  ~PetscRunner()
70  {
71  exit();
72  }
73 
75  void init(M const& A) override
76  {
77  exit();
78 #if HAVE_MPI
79  KSPCreate(Environment::comm(), &ksp_);
80 #else
81  KSPCreate(PETSC_COMM_SELF, &ksp_);
82 #endif
83 
84  PETSc::KSPSetOperators(ksp_, A, A);
85  initKSP(ksp_, prefix_);
86  initialized_ = true;
87  }
88 
90  void exit() override
91  {
92  if (initialized_)
93  KSPDestroy(&ksp_);
94  initialized_ = false;
95  }
96 
98  int solve([[maybe_unused]] M const& A, X& x, Y const& b, SolverInfo& solverInfo) override
99  {
100  assert(initialized_);
101 
102  KSPSolve(ksp_, b, x);
103 
104  if (info_ > 0)
105  KSPReasonView(ksp_, PETSC_VIEWER_STDOUT_WORLD);
106 
107  KSPConvergedReason reason;
108  KSPGetConvergedReason(ksp_, &reason);
109  return reason > 0 ? 0 : int(reason);
110  }
111 
112  KSP ksp() { return ksp_; }
113 
114  protected:
115  // initialize the KSP solver from the initfile
116  virtual void initKSP(KSP ksp, std::string prefix) const
117  {
118  // see https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPType.html
119  auto kspType = Parameters::get<std::string>(prefix + "->ksp");
120  std::string kspTypeStr = kspType.value_or("default");
121 
122  if (!kspType)
123  Parameters::get(prefix, kspTypeStr);
124 
125  if (kspTypeStr == "direct")
126  initDirectSolver(ksp, prefix);
127  else if (kspTypeStr != "default") {
128  KSPSetType(ksp, kspTypeStr.c_str());
129 
130  // initialize some KSP specific parameters
131  initKSPParameters(ksp, kspTypeStr.c_str(), prefix);
132 
133  // set initial guess to nonzero only for non-preonly ksp type
134  if (kspTypeStr != "preonly")
135  KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
136  }
137 
138  // set a KSPMonitor if info > 0
139  int info = 0;
140  Parameters::get(prefix + "->info", info);
141  if (info == 1) {
142  KSPMonitorSet(ksp, PETSc::KSPMonitorCyclic, nullptr, nullptr);
143  } else if (info > 1) {
144  KSPMonitorSet(ksp, PETSc::KSPMonitorNoisy, nullptr, nullptr);
145  }
146 
147  // set solver tolerances
148  auto maxIt = Parameters::get<PetscInt>(prefix + "->max it");
149  auto rTol = Parameters::get<PetscReal>(prefix + "->rtol");
150  auto aTol = Parameters::get<PetscReal>(prefix + "->atol");
151  auto dTol = Parameters::get<PetscReal>(prefix + "->divtol");
152  if (maxIt || rTol || aTol || dTol)
153  KSPSetTolerances(ksp,
154  rTol.value_or(PETSC_DEFAULT), aTol.value_or(PETSC_DEFAULT),
155  dTol.value_or(PETSC_DEFAULT), maxIt.value_or(PETSC_DEFAULT));
156 
157  auto kspPrefixParam = Parameters::get<std::string>(prefix + "->prefix");
158  if (kspPrefixParam) {
159  std::string kspPrefix = kspPrefixParam.value() + "_";
160  KSPSetOptionsPrefix(ksp, kspPrefix.c_str());
161  KSPSetFromOptions(ksp);
162  }
163 
164  PC pc;
165  KSPGetPC(ksp, &pc);
166  initPC(pc, prefix + "->pc");
167 
168  // show details of ksp object
169  if (info >= 10) {
170  KSPView(ksp, PETSC_VIEWER_STDOUT_WORLD);
171  }
172  }
173 
174 
175  // initialize a direct solver from the initfile
176  virtual void initDirectSolver(KSP ksp, std::string prefix) const
177  {
178  KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
179  KSPSetType(ksp, KSPRICHARDSON);
180 
181  PC pc;
182  KSPGetPC(ksp, &pc);
183  PCSetType(pc, PCLU);
184  auto matSolverType = Parameters::get<std::string>(prefix + "->mat solver");
185  if (matSolverType)
186  PETSc::PCFactorSetMatSolverType(pc, matSolverType.value().c_str());
187  else {
188  Mat Amat, Pmat;
189  KSPGetOperators(ksp, &Amat, &Pmat);
190  MatType type;
191  MatGetType(Pmat, &type);
192  if (std::strcmp(type, MATSEQAIJ) == 0)
193  PETSc::PCFactorSetMatSolverType(pc, MATSOLVERUMFPACK);
194  else if (std::strcmp(type, MATAIJ) == 0)
195  PETSc::PCFactorSetMatSolverType(pc, MATSOLVERMUMPS);
196  }
197  }
198 
199 
200  // initialize the preconditioner pc from the initfile
201  virtual void initPC(PC pc, std::string prefix) const
202  {
203  // see https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCType.html
204  auto pcType = Parameters::get<std::string>(prefix);
205  std::string pcTypeStr = "default";
206  if (pcType) {
207  pcTypeStr = pcType.value();
208  PCSetType(pc, pcTypeStr.c_str());
209  }
210 
211  if (pcTypeStr == "lu") {
212  // configure matsolverpackage
213  // see https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSolverType.html
214  auto matSolverType = Parameters::get<std::string>(prefix + "->mat solver");
215  if (matSolverType)
216  PETSc::PCFactorSetMatSolverType(pc, matSolverType.value().c_str());
217  } else if (pcTypeStr == "bjacobi") {
218  // configure sub solver
219  PetscInt n_local, first;
220  KSP* ksps;
221  PCSetUp(pc);
222  PCBJacobiGetSubKSP(pc, &n_local, &first, &ksps);
223  for (PetscInt i = 0; i < n_local; ++i)
224  initKSP(ksps[i], prefix + "->sub ksp");
225  } else if (pcTypeStr == "ksp") {
226  // configure sub ksp type
227  KSP ksp;
228  PCSetUp(pc);
229  PCKSPGetKSP(pc, &ksp);
230  initKSP(ksp, prefix + "->ksp");
231  }
232 
233  auto pcPrefixParam = Parameters::get<std::string>(prefix + "->prefix");
234  if (pcPrefixParam) {
235  std::string pcPrefix = pcPrefixParam.value() + "_";
236  PCSetOptionsPrefix(pc, pcPrefix.c_str());
237  PCSetFromOptions(pc);
238  }
239  }
240 
241 
242  // provide initfile parameters for some PETSc KSP parameters
243  virtual void initKSPParameters(KSP ksp, char const* ksptype, std::string prefix) const
244  {
245  // parameters for the Richardson solver
246  if (std::strcmp(ksptype, KSPRICHARDSON) == 0)
247  {
248  auto scale = Parameters::get<PetscReal>(prefix + "->scale");
249  if (scale)
250  KSPRichardsonSetScale(ksp, scale.value());
251  }
252  // parameters for the gmres solver
253  else if (std::strcmp(ksptype, KSPGMRES) == 0 ||
254  std::strcmp(ksptype, KSPFGMRES) == 0 ||
255  std::strcmp(ksptype, KSPLGMRES) == 0 ||
256  std::strcmp(ksptype, KSPDGGMRES) == 0 ||
257  std::strcmp(ksptype, KSPPGMRES) == 0)
258  {
259  auto restart = Parameters::get<PetscInt>(prefix + "->restart");
260  if (restart)
261  KSPGMRESSetRestart(ksp, restart.value());
262 
263  auto gramschmidt = Parameters::get<std::string>(prefix + "->gramschmidt");
264  if (gramschmidt) {
265  if (gramschmidt.value() == "classical") {
266  KSPGMRESSetOrthogonalization(ksp, KSPGMRESClassicalGramSchmidtOrthogonalization);
267  auto refinementType = Parameters::get<std::string>(prefix + "->refinement type");
268  if (refinementType) {
269  if (refinementType.value() == "never")
270  KSPGMRESSetCGSRefinementType(ksp, KSP_GMRES_CGS_REFINE_NEVER);
271  else if (refinementType.value() == "ifneeded")
272  KSPGMRESSetCGSRefinementType(ksp, KSP_GMRES_CGS_REFINE_IFNEEDED);
273  else if (refinementType.value() == "always")
274  KSPGMRESSetCGSRefinementType(ksp, KSP_GMRES_CGS_REFINE_ALWAYS);
275  }
276  } else if (gramschmidt.value() == "modified") {
277  KSPGMRESSetOrthogonalization(ksp, KSPGMRESModifiedGramSchmidtOrthogonalization);
278  }
279  }
280  }
281  // parameters for the BiCGStab_ell solver
282  else if (std::strcmp(ksptype, KSPBCGSL) == 0)
283  {
284  auto ell = Parameters::get<PetscInt>(prefix + "->ell");
285  if (ell)
286  KSPBCGSLSetEll(ksp, ell.value());
287  }
288  // parameters for the GCR solver
289  else if (std::strcmp(ksptype, KSPGCR) == 0)
290  {
291  auto restart = Parameters::get<PetscInt>(prefix + "->restart");
292  if (restart)
293  KSPGCRSetRestart(ksp, restart.value());
294  }
295  }
296 
297  protected:
298  std::string prefix_;
299  int info_ = 0;
300 
301  KSP ksp_;
302  bool initialized_ = false;
303  };
304 
305 } // end namespace AMDiS
PetscRunner(std::string const &prefix)
Constructor.
Definition: PetscRunner.hpp:63
Wrapper around PETSc KSP and PC objects to solve a linear system.
Definition: PetscRunner.hpp:46
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
void init(M const &A) override
Implements RunnerInterface::init()
Definition: PetscRunner.hpp:75
int solve([[maybe_unused]] M const &A, X &x, Y const &b, SolverInfo &solverInfo) override
Implements RunnerInterface::solve()
Definition: PetscRunner.hpp:98
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
void info(int noInfoLevel, std::string const &str, Args &&... args)
prints a message, if Environment::infoLevel() >= noInfoLevel
Definition: Output.hpp:105
Definition: SolverInfo.hpp:11
void exit() override
Implements RunnerInterface::exit()
Definition: PetscRunner.hpp:90
Interface for Runner / Worker types used in solver classes.
Definition: RunnerInterface.hpp:11
static Dune::MPIHelper::MPICommunicator comm()
Return the MPI_Comm object (or a fake communicator)
Definition: Environment.hpp:86