/usr/include/dolfin/nls/TAOLinearBoundSolver.h is in libdolfin-dev 1.4.0+dfsg-4.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 | // Copyright (C) 2012 Corrado Maurini
//
// This file is part of DOLFIN.
//
// DOLFIN is free software: you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// DOLFIN is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public License
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
//
// First added: 2012-12-03
// Last changed: 2013-04-02
#ifndef _TAOLinearBoundSolver_H
#define _TAOLinearBoundSolver_H
#ifdef ENABLE_PETSC_TAO
#include <memory>
#include <dolfin/log/dolfin_log.h>
#include <dolfin/common/NoDeleter.h>
#include <dolfin/common/types.h>
#include <map>
#include <petscksp.h>
#include <petscpc.h>
#include <petsctao.h>
#include <dolfin/la/PETScObject.h>
#include <dolfin/la/KrylovSolver.h>
namespace dolfin
{
/// Forward declarations
class GenericMatrix;
class GenericVector;
class PETScMatrix;
class PETScVector;
class PETScPreconditioner;
class PETScKrylovSolver;
class PETScKSPDeleter;
/// This class provides a bound constrained solver for a
/// linear variational inequality defined by a matrix A and a vector b.
/// It solves the problem:
///
/// Find :math:`x_l\leq x\leq x_u` such that
/// :math:`(Ax-b)\cdot (y-x)\geq 0,\; \forall x_l\leq y\leq x_u`
///
/// It is a wrapper for the TAO bound constrained solver.
///
/// *Example*
/// .. code-block:: python
///
/// # Assemble the linear system
/// A, b = assemble_system(a, L, bc)
/// # Define the constraints
/// constraint_u = Constant(1.)
/// constraint_l = Constant(0.)
/// u_min = interpolate(constraint_l, V)
/// u_max = interpolate(constraint_u, V)
/// # Define the function to store the solution
/// usol=Function(V)
/// # Create the TAOLinearBoundSolver
/// solver=TAOLinearBoundSolver("tao_gpcg","gmres")
/// # Set some parameters
/// solver.parameters["monitor_convergence"]=True
/// solver.parameters["report"]=True
/// # Solve the problem
/// solver.solve(A, usol.vector(), b , u_min.vector(), u_max.vector())
/// info(solver.parameters,True)
///
class TAOLinearBoundSolver : public Variable, public PETScObject
{
public:
/// Create TAO bound constrained solver
TAOLinearBoundSolver(const std::string method = "default",
const std::string ksp_type = "default",
const std::string pc_type = "default");
/// Destructor
~TAOLinearBoundSolver();
/// Solve the linear variational inequality defined by A and b
/// with xl =< x <= xu
std::size_t solve(const GenericMatrix& A, GenericVector& x,
const GenericVector& b, const GenericVector& xl,
const GenericVector& xu);
/// Solve the linear variational inequality defined by A and b
/// with xl =< x <= xu
std::size_t solve(const PETScMatrix& A, PETScVector& x,
const PETScVector& b,
const PETScVector& xl, const PETScVector& xu);
// Set the TAO solver type
void set_solver(const std::string&);
/// Set PETSC Krylov Solver (ksp) used by TAO
void set_ksp( const std::string ksp_type = "default");
// Return TAO solver pointer
Tao tao() const;
/// Return a list of available Tao solver methods
static std::vector<std::pair<std::string, std::string> > methods();
/// Return a list of available krylov solvers
static std::vector<std::pair<std::string, std::string> > krylov_solvers();
/// Return a list of available preconditioners
static std::vector<std::pair<std::string, std::string> > preconditioners();
/// Default parameter values
static Parameters default_parameters()
{
Parameters p("tao_solver");
p.add("monitor_convergence" , false);
p.add("report" , false);
p.add("function_absolute_tol" , 1.0e-10);
p.add("function_relative_tol" , 1.0e-10);
p.add("gradient_absolute_tol" , 1.0e-8);
p.add("gradient_relative_tol" , 1.0e-8);
p.add("gradient_t_tol" , 0.0);
p.add("error_on_nonconvergence", true);
p.add("maximum_iterations" , 100);
p.add("options_prefix" , "default");
Parameters ksp("krylov_solver");
ksp = KrylovSolver::default_parameters();
p.add(ksp);
return p;
}
// Return Matrix shared pointer
std::shared_ptr<const PETScMatrix> get_matrix() const;
// Return load vector shared pointer
std::shared_ptr<const PETScVector> get_vector() const;
private:
// Set operators with GenericMatrix and GenericVector
void set_operators(std::shared_ptr<const GenericMatrix> A,
std::shared_ptr<const GenericVector> b);
// Set operators with shared pointer to PETSc objects
void set_operators(std::shared_ptr<const PETScMatrix> A,
std::shared_ptr<const PETScVector> b);
// Callback for changes in parameter values
void read_parameters();
// Available ksp solvers
static const std::map<std::string, const KSPType> _ksp_methods;
// Available tao solvers descriptions
static const std::vector<std::pair<std::string, std::string> >
_methods_descr;
// Set options
void set_ksp_options();
// Initialize TAO solver
void init(const std::string& method);
// Tao solver pointer
Tao _tao;
// Petsc preconditioner
std::shared_ptr<PETScPreconditioner> preconditioner;
// Operator (the matrix) and the vector
std::shared_ptr<const PETScMatrix> A;
std::shared_ptr<const PETScVector> b;
bool preconditioner_set;
// Computes the value of the objective function and its gradient.
static PetscErrorCode
__TAOFormFunctionGradientQuadraticProblem(Tao tao, Vec X,
PetscReal *ener, Vec G,
void *ptr);
// Computes the hessian of the quadratic objective function
static PetscErrorCode
__TAOFormHessianQuadraticProblem(Tao tao,Vec X, Mat H, Mat Hpre,
void *ptr);
//-------------------------------------------------------------------------
// Monitor the state of the solution at each iteration. The
// output printed to the screen is:
//
// iterate - the current iterate number (>=0)
// f - the current function value
// gnorm - the square of the gradient norm, duality gap, or other
// measure
// indicating distance from optimality.
// cnorm - the infeasibility of the current solution with regard
// to the constraints.
// xdiff - the step length or trust region radius of the most
// recent iterate.
//-------------------------------------------------------------------------
static PetscErrorCode __TAOMonitor(Tao tao, void *ctx);
};
}
#endif
#endif
|