Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
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/solver_control.h>
26 # include <deal.II/lac/vector.h>
27 # include <deal.II/lac/la_parallel_vector.h>
28 
29 # include <memory>
30 
31 # include <Epetra_LinearProblem.h>
32 # include <AztecOO.h>
33 # include <Epetra_Operator.h>
34 # include <Amesos.h>
35 
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 namespace TrilinosWrappers
40 {
41  // forward declarations
42  class SparseMatrix;
43  class PreconditionBase;
44 
45 
65  class SolverBase
66  {
67  public:
68 
77  {
81  cg,
85  cgs,
98  } solver_name;
99 
105  {
114  explicit
115  AdditionalData (const bool output_solver_details = false,
116  const unsigned int gmres_restart_parameter = 30);
117 
123 
127  const unsigned int gmres_restart_parameter;
128  };
129 
134  const AdditionalData &data = AdditionalData());
135 
140  SolverBase (const enum SolverName solver_name,
141  SolverControl &cn,
142  const AdditionalData &data = AdditionalData());
143 
147  virtual ~SolverBase () = default;
148 
154  void
155  solve (const SparseMatrix &A,
156  MPI::Vector &x,
157  const MPI::Vector &b,
158  const PreconditionBase &preconditioner);
159 
167  void
168  solve (const Epetra_Operator &A,
169  MPI::Vector &x,
170  const MPI::Vector &b,
171  const PreconditionBase &preconditioner);
172 
182  void
183  solve (const Epetra_Operator &A,
184  MPI::Vector &x,
185  const MPI::Vector &b,
186  const Epetra_Operator &preconditioner);
187 
197  void
198  solve (const Epetra_Operator &A,
199  Epetra_MultiVector &x,
200  const Epetra_MultiVector &b,
201  const PreconditionBase &preconditioner);
202 
213  void
214  solve (const Epetra_Operator &A,
215  Epetra_MultiVector &x,
216  const Epetra_MultiVector &b,
217  const Epetra_Operator &preconditioner);
218 
219 
220 
231  void
232  solve (const SparseMatrix &A,
233  ::Vector<double> &x,
234  const ::Vector<double> &b,
235  const PreconditionBase &preconditioner);
236 
248  void
249  solve (Epetra_Operator &A,
250  ::Vector<double> &x,
251  const ::Vector<double> &b,
252  const PreconditionBase &preconditioner);
253 
260  void
261  solve (const SparseMatrix &A,
263  const ::LinearAlgebra::distributed::Vector<double> &b,
264  const PreconditionBase &preconditioner);
265 
273  void
274  solve (Epetra_Operator &A,
276  const ::LinearAlgebra::distributed::Vector<double> &b,
277  const PreconditionBase &preconditioner);
278 
279 
283  SolverControl &control() const;
284 
289  int,
290  << "An error with error number " << arg1
291  << " occurred while calling a Trilinos function");
292 
293  protected:
294 
302 
303  private:
304 
309  template <typename Preconditioner>
310  void do_solve(const Preconditioner &preconditioner);
311 
315  template <typename Preconditioner>
316  void set_preconditioner (AztecOO &solver,
317  const Preconditioner &preconditioner);
318 
324  std::unique_ptr<Epetra_LinearProblem> linear_problem;
325 
330  std::unique_ptr<AztecOO_StatusTest> status_test;
331 
336  AztecOO solver;
337 
342  };
343 
344 
345  // provide a declaration for two explicit specializations
346  template <>
347  void
348  SolverBase::set_preconditioner(AztecOO &solver,
349  const PreconditionBase &preconditioner);
350 
351  template <>
352  void
353  SolverBase::set_preconditioner(AztecOO &solver,
354  const Epetra_Operator &preconditioner);
355 
356 
357 
364  class SolverCG : public SolverBase
365  {
366  public:
372  {
376  explicit
377  AdditionalData (const bool output_solver_details = false);
378  };
379 
387  SolverCG (SolverControl &cn,
388  const AdditionalData &data = AdditionalData());
389 
390  protected:
395  };
396 
397 
398 
405  class SolverCGS : public SolverBase
406  {
407  public:
412  {
416  explicit
417  AdditionalData (const bool output_solver_details = false);
418  };
419 
428  const AdditionalData &data = AdditionalData());
429 
430  protected:
435  };
436 
437 
438 
445  class SolverGMRES : public SolverBase
446  {
447  public:
452  {
457  explicit
458  AdditionalData (const bool output_solver_details = false,
459  const unsigned int restart_parameter = 30);
460  };
461 
470  const AdditionalData &data = AdditionalData());
471 
472  protected:
477  };
478 
479 
480 
488  class SolverBicgstab : public SolverBase
489  {
490  public:
495  {
499  explicit
500  AdditionalData (const bool output_solver_details = false);
501  };
502 
511  const AdditionalData &data = AdditionalData());
512 
513  protected:
518  };
519 
520 
521 
529  class SolverTFQMR : public SolverBase
530  {
531  public:
536  {
540  explicit
541  AdditionalData (const bool output_solver_details = false);
542  };
543 
552  const AdditionalData &data = AdditionalData());
553 
554  protected:
559  };
560 
561 
562 
577  {
578  public:
579 
585  {
589  explicit
590  AdditionalData (const bool output_solver_details = false,
591  const std::string &solver_type = "Amesos_Klu");
592 
598 
617  std::string solver_type;
618  };
619 
624  const AdditionalData &data = AdditionalData());
625 
629  virtual ~SolverDirect () = default;
630 
637  void initialize (const SparseMatrix &A);
638 
644  void solve (MPI::Vector &x, const MPI::Vector &b);
645 
651  void
653  const ::LinearAlgebra::distributed::Vector<double> &b);
654 
661  void
662  solve (const SparseMatrix &A,
663  MPI::Vector &x,
664  const MPI::Vector &b);
665 
673  void
674  solve (const SparseMatrix &A,
675  ::Vector<double> &x,
676  const ::Vector<double> &b);
677 
684  void
685  solve (const SparseMatrix &A,
687  const ::LinearAlgebra::distributed::Vector<double> &b);
688 
692  SolverControl &control() const;
693 
698  int,
699  << "An error with error number " << arg1
700  << " occurred while calling a Trilinos function");
701 
702  private:
707  void do_solve();
708 
716 
722  std::unique_ptr<Epetra_LinearProblem> linear_problem;
723 
728  std::unique_ptr<Amesos_BaseSolver> solver;
729 
734  };
735 
736 }
737 
738 DEAL_II_NAMESPACE_CLOSE
739 
740 #endif // DEAL_II_WITH_TRILINOS
741 
742 /*---------------------------- trilinos_solver.h ---------------------------*/
743 
744 #endif
745 /*---------------------------- trilinos_solver.h ---------------------------*/
void set_preconditioner(AztecOO &solver, const Preconditioner &preconditioner)
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
virtual ~SolverDirect()=default
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)
const AdditionalData additional_data
void do_solve(const Preconditioner &preconditioner)
virtual ~SolverBase()=default
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:346
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
void solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner)
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverTFQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
AdditionalData(const bool output_solver_details=false)
SolverBase(SolverControl &cn, const AdditionalData &data=AdditionalData())
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
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)