18 #ifndef _TAOLinearBoundSolver_H 19 #define _TAOLinearBoundSolver_H 29 #include <dolfin/common/NoDeleter.h> 30 #include <dolfin/common/types.h> 32 #include <dolfin/la/PETScObject.h> 33 #include <dolfin/la/KrylovSolver.h> 43 class PETScPreconditioner;
44 class PETScKrylovSolver;
45 class PETScKSPDeleter;
85 const std::string ksp_type =
"default",
86 const std::string pc_type =
"default");
107 void set_ksp(
const std::string ksp_type =
"default");
113 static std::map<std::string, std::string>
methods();
126 p.
add(
"monitor_convergence" ,
false);
127 p.
add(
"report" ,
false);
128 p.
add(
"gradient_absolute_tol" , 1.0e-8);
129 p.
add(
"gradient_relative_tol" , 1.0e-8);
130 p.
add(
"gradient_t_tol" , 0.0);
131 p.
add(
"error_on_nonconvergence",
true);
132 p.
add(
"maximum_iterations" , 100);
133 p.
add(
"options_prefix" ,
"default");
143 std::shared_ptr<const PETScMatrix>
get_matrix()
const;
146 std::shared_ptr<const PETScVector>
get_vector()
const;
151 void set_operators(std::shared_ptr<const GenericMatrix> A,
152 std::shared_ptr<const GenericVector> b);
155 void set_operators(std::shared_ptr<const PETScMatrix> A,
156 std::shared_ptr<const PETScVector> b);
159 void read_parameters();
162 static const std::map<std::string, const KSPType> _ksp_methods;
165 static const std::map<std::string, std::string> _methods_descr;
168 void set_ksp_options();
171 void init(
const std::string& method);
177 std::shared_ptr<PETScPreconditioner> _preconditioner;
180 std::shared_ptr<const PETScMatrix> _matA;
181 std::shared_ptr<const PETScVector> _b;
183 bool _preconditioner_set;
186 static PetscErrorCode
187 __TAOFormFunctionGradientQuadraticProblem(Tao tao, Vec X,
188 PetscReal *ener, Vec G,
192 static PetscErrorCode
193 __TAOFormHessianQuadraticProblem(Tao tao,Vec X, Mat H, Mat Hpre,
210 static PetscErrorCode __TAOMonitor(Tao tao,
void *ctx);
void set_solver(const std::string &)
Set the TAO solver type.
Definition: TAOLinearBoundSolver.cpp:272
Tao tao() const
Return TAO solver pointer.
Definition: TAOLinearBoundSolver.cpp:327
Definition: PETScVector.h:60
Common base class for DOLFIN variables.
Definition: Variable.h:35
std::shared_ptr< const PETScVector > get_vector() const
Return load vector shared pointer.
Definition: TAOLinearBoundSolver.cpp:337
static std::map< std::string, std::string > preconditioners()
Return a list of available preconditioners.
Definition: TAOLinearBoundSolver.cpp:74
void add(std::string key)
Definition: Parameters.h:128
Definition: TAOLinearBoundSolver.h:76
~TAOLinearBoundSolver()
Destructor.
Definition: TAOLinearBoundSolver.cpp:119
Definition: PETScMatrix.h:58
TAOLinearBoundSolver(MPI_Comm comm)
Create TAO bound constrained solver.
Definition: TAOLinearBoundSolver.cpp:79
std::shared_ptr< const PETScMatrix > get_matrix() const
Return Matrix shared pointer.
Definition: TAOLinearBoundSolver.cpp:332
static std::map< std::string, std::string > methods()
Return a list of available Tao solver methods.
Definition: TAOLinearBoundSolver.cpp:64
Definition: Parameters.h:94
std::size_t solve(const GenericMatrix &A, GenericVector &x, const GenericVector &b, const GenericVector &xl, const GenericVector &xu)
Definition: TAOLinearBoundSolver.cpp:142
This class defines a common interface for matrices.
Definition: GenericMatrix.h:46
static std::map< std::string, std::string > krylov_solvers()
Return a list of available krylov solvers.
Definition: TAOLinearBoundSolver.cpp:69
void set_ksp(const std::string ksp_type="default")
Set PETSC Krylov Solver (ksp) used by TAO.
Definition: TAOLinearBoundSolver.cpp:302
Definition: PETScObject.h:33
static Parameters default_parameters()
Default parameter values.
Definition: KrylovSolver.cpp:32
This class defines a common interface for vectors.
Definition: GenericVector.h:47
static Parameters default_parameters()
Default parameter values.
Definition: TAOLinearBoundSolver.h:122