Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
trilinos_solver.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_trilinos_solver_h
17 # define dealii_trilinos_solver_h
18 
19 
20 # include <deal.II/base/config.h>
21 
22 # ifdef DEAL_II_WITH_TRILINOS
23 
24 # include <deal.II/lac/exceptions.h>
25 # include <deal.II/lac/la_parallel_vector.h>
26 # include <deal.II/lac/solver_control.h>
27 # include <deal.II/lac/vector.h>
28 
29 # include <Amesos.h>
30 # include <AztecOO.h>
31 # include <Epetra_LinearProblem.h>
32 # include <Epetra_Operator.h>
33 
34 # include <memory>
35 
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 namespace TrilinosWrappers
40 {
41  // forward declarations
42  class SparseMatrix;
43  class PreconditionBase;
44 
45 
66  class SolverBase
67  {
68  public:
77  {
81  cg,
85  cgs,
98  } solver_name;
99 
105  {
114  explicit AdditionalData(const bool output_solver_details = false,
115  const unsigned int gmres_restart_parameter = 30);
116 
122 
126  const unsigned int gmres_restart_parameter;
127  };
128 
133  const AdditionalData &data = AdditionalData());
134 
139  SolverBase(const enum SolverName solver_name,
140  SolverControl & cn,
141  const AdditionalData &data = AdditionalData());
142 
146  virtual ~SolverBase() = default;
147 
153  void
154  solve(const SparseMatrix & A,
155  MPI::Vector & x,
156  const MPI::Vector & b,
157  const PreconditionBase &preconditioner);
158 
166  void
167  solve(const Epetra_Operator & A,
168  MPI::Vector & x,
169  const MPI::Vector & b,
170  const PreconditionBase &preconditioner);
171 
181  void
182  solve(const Epetra_Operator &A,
183  MPI::Vector & x,
184  const MPI::Vector & b,
185  const Epetra_Operator &preconditioner);
186 
196  void
197  solve(const Epetra_Operator & A,
198  Epetra_MultiVector & x,
199  const Epetra_MultiVector &b,
200  const PreconditionBase & preconditioner);
201 
212  void
213  solve(const Epetra_Operator & A,
214  Epetra_MultiVector & x,
215  const Epetra_MultiVector &b,
216  const Epetra_Operator & preconditioner);
217 
218 
219 
230  void
231  solve(const SparseMatrix & A,
232  ::Vector<double> & x,
233  const ::Vector<double> &b,
234  const PreconditionBase & preconditioner);
235 
247  void
248  solve(Epetra_Operator & A,
249  ::Vector<double> & x,
250  const ::Vector<double> &b,
251  const PreconditionBase & preconditioner);
252 
259  void
260  solve(const SparseMatrix & A,
262  const ::LinearAlgebra::distributed::Vector<double> &b,
263  const PreconditionBase &preconditioner);
264 
272  void
273  solve(Epetra_Operator & A,
275  const ::LinearAlgebra::distributed::Vector<double> &b,
276  const PreconditionBase &preconditioner);
277 
278 
282  SolverControl &
283  control() const;
284 
289  int,
290  << "An error with error number " << arg1
291  << " occurred while calling a Trilinos function");
292 
293  protected:
301 
302  private:
307  template <typename Preconditioner>
308  void
309  do_solve(const Preconditioner &preconditioner);
310 
314  template <typename Preconditioner>
315  void
316  set_preconditioner(AztecOO &solver, const Preconditioner &preconditioner);
317 
323  std::unique_ptr<Epetra_LinearProblem> linear_problem;
324 
329  std::unique_ptr<AztecOO_StatusTest> status_test;
330 
335  AztecOO solver;
336 
341  };
342 
343 
344  // provide a declaration for two explicit specializations
345  template <>
346  void
347  SolverBase::set_preconditioner(AztecOO & solver,
348  const PreconditionBase &preconditioner);
349 
350  template <>
351  void
352  SolverBase::set_preconditioner(AztecOO & solver,
353  const Epetra_Operator &preconditioner);
354 
355 
356 
363  class SolverCG : public SolverBase
364  {
365  public:
371  {
375  explicit AdditionalData(const bool output_solver_details = false);
376  };
377 
385  SolverCG(SolverControl &cn, const AdditionalData &data = AdditionalData());
386 
387  protected:
392  };
393 
394 
395 
402  class SolverCGS : public SolverBase
403  {
404  public:
409  {
413  explicit AdditionalData(const bool output_solver_details = false);
414  };
415 
424 
425  protected:
430  };
431 
432 
433 
440  class SolverGMRES : public SolverBase
441  {
442  public:
447  {
452  explicit AdditionalData(const bool output_solver_details = false,
453  const unsigned int restart_parameter = 30);
454  };
455 
464  const AdditionalData &data = AdditionalData());
465 
466  protected:
471  };
472 
473 
474 
482  class SolverBicgstab : public SolverBase
483  {
484  public:
489  {
493  explicit AdditionalData(const bool output_solver_details = false);
494  };
495 
504  const AdditionalData &data = AdditionalData());
505 
506  protected:
511  };
512 
513 
514 
522  class SolverTFQMR : public SolverBase
523  {
524  public:
529  {
533  explicit AdditionalData(const bool output_solver_details = false);
534  };
535 
544  const AdditionalData &data = AdditionalData());
545 
546  protected:
551  };
552 
553 
554 
569  {
570  public:
576  {
580  explicit AdditionalData(const bool output_solver_details = false,
581  const std::string &solver_type = "Amesos_Klu");
582 
588 
607  std::string solver_type;
608  };
609 
614  const AdditionalData &data = AdditionalData());
615 
619  virtual ~SolverDirect() = default;
620 
627  void
628  initialize(const SparseMatrix &A);
629 
635  void
636  solve(MPI::Vector &x, const MPI::Vector &b);
637 
643  void
645  const ::LinearAlgebra::distributed::Vector<double> &b);
646 
653  void
654  solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b);
655 
663  void
664  solve(const SparseMatrix & A,
665  ::Vector<double> & x,
666  const ::Vector<double> &b);
667 
674  void
675  solve(const SparseMatrix & A,
677  const ::LinearAlgebra::distributed::Vector<double> &b);
678 
682  SolverControl &
683  control() const;
684 
689  int,
690  << "An error with error number " << arg1
691  << " occurred while calling a Trilinos function");
692 
693  private:
698  void
699  do_solve();
700 
708 
714  std::unique_ptr<Epetra_LinearProblem> linear_problem;
715 
720  std::unique_ptr<Amesos_BaseSolver> solver;
721 
726  };
727 
728 } // namespace TrilinosWrappers
729 
730 DEAL_II_NAMESPACE_CLOSE
731 
732 # endif // DEAL_II_WITH_TRILINOS
733 
734 /*---------------------------- trilinos_solver.h ---------------------------*/
735 
736 #endif
737 /*---------------------------- trilinos_solver.h ---------------------------*/
void set_preconditioner(AztecOO &solver, const Preconditioner &preconditioner)
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
virtual ~SolverDirect()=default
SolverBase(SolverControl &cn, const AdditionalData &data=AdditionalData())
std::unique_ptr< Amesos_BaseSolver > solver
const AdditionalData additional_data
std::unique_ptr< AztecOO_StatusTest > status_test
const AdditionalData additional_data
AdditionalData(const bool output_solver_details=false, const unsigned int restart_parameter=30)
const AdditionalData additional_data
void initialize(const SparseMatrix &A)
static ::ExceptionBase & ExcTrilinosError(int arg1)
void solve(MPI::Vector &x, const MPI::Vector &b)
void solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner)
const AdditionalData additional_data
virtual ~SolverBase()=default
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:518
SolverCGS(SolverControl &cn, const AdditionalData &data=AdditionalData())
AdditionalData(const bool output_solver_details=false)
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverControl & control() const
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverTFQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
AdditionalData(const bool output_solver_details=false)
const AdditionalData additional_data
AdditionalData(const bool output_solver_details=false, const std::string &solver_type="Amesos_Klu")
SolverDirect(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
AdditionalData(const bool output_solver_details=false)
SolverControl & control() const
void do_solve(const Preconditioner &preconditioner)
AdditionalData(const bool output_solver_details=false, const unsigned int gmres_restart_parameter=30)
SolverGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
static ::ExceptionBase & ExcTrilinosError(int arg1)
AdditionalData(const bool output_solver_details=false)