dolfin/nls

Documentation for C++ code found in dolfin/nls/*.h

Classes

NewtonSolver

C++ documentation for NewtonSolver from dolfin/nls/NewtonSolver.h:

class dolfin::NewtonSolver

This class defines a Newton solver for nonlinear systems of equations of the form \(F(x) = 0\) .

dolfin::NewtonSolver::NewtonSolver(MPI_Comm comm, std::shared_ptr<GenericLinearSolver> solver, GenericLinearAlgebraFactory &factory)

Create nonlinear solver using provided linear solver Arguments comm (MPI_Ccmm) The MPI communicator. solver (GenericLinearSolver ) The linear solver. factory (GenericLinearAlgebraFactory ) The factory.

Parameters:
  • comm
  • solver
  • factory
dolfin::NewtonSolver::NewtonSolver(MPI_Comm comm = MPI_COMM_WORLD)

Create nonlinear solver.

Parameters:comm
bool dolfin::NewtonSolver::converged(const GenericVector &r, const NonlinearProblem &nonlinear_problem, std::size_t iteration)

Convergence test. It may be overloaded using virtual inheritance and this base criterion may be called from derived, both in C++ and Python. Arguments r (GenericVector ) Residual for criterion evaluation. nonlinear_problem (NonlinearProblem ) The nonlinear problem. iteration (std::size_t) Newton iteration number. Returns bool Whether convergence occurred.

Parameters:
  • r
  • nonlinear_problem
  • iteration
double dolfin::NewtonSolver::get_relaxation_parameter()

Get relaxation parameter Returns double Relaxation parameter value.

std::size_t dolfin::NewtonSolver::iteration() const

Return current Newton iteration number Returns std::size_t The iteration number.

std::size_t dolfin::NewtonSolver::krylov_iterations() const

Return number of Krylov iterations elapsed since solve started Returns std::size_t The number of iterations.

GenericLinearSolver &dolfin::NewtonSolver::linear_solver() const

Return the linear solver Returns:cpp:any:GenericLinearSolver The linear solver.

double dolfin::NewtonSolver::relative_residual() const

Return current relative residual Returns double Current relative residual.

double dolfin::NewtonSolver::residual() const

Return current residual Returns double Current residual.

double dolfin::NewtonSolver::residual0() const

Return initial residual Returns double Initial residual.

void dolfin::NewtonSolver::set_relaxation_parameter(double relaxation_parameter)

Set relaxation parameter. Default value 1.0 means full Newton method, value smaller than 1.0 relaxes the method by shrinking effective Newton step size by the given factor. Arguments relaxation_parameter(double) Relaxation parameter value.

Parameters:relaxation_parameter
std::pair<std::size_t, bool> dolfin::NewtonSolver::solve(NonlinearProblem &nonlinear_function, GenericVector &x)

Solve abstract nonlinear problem \(F(x) = 0\) for given \(F\) and Jacobian \(\dfrac{\partial F}{\partial x}\) . Arguments nonlinear_function (NonlinearProblem ) The nonlinear problem. x (GenericVector ) The vector. Returns std::pair<std::size_t, bool> Pair of number of Newton iterations, and whether iteration converged)

Parameters:
  • nonlinear_function
  • x
void dolfin::NewtonSolver::solver_setup(std::shared_ptr<const GenericMatrix> A, std::shared_ptr<const GenericMatrix> P, const NonlinearProblem &nonlinear_problem, std::size_t iteration)

Setup solver to be used with system matrix A and preconditioner matrix P. It may be overloaded to get finer control over linear solver setup, various linesearch tricks, etc. Note that minimal implementation should call set_operators method of the linear solver. Arguments A (std::shared_ptr<const GenericMatrix>) System Jacobian matrix. J (std::shared_ptr<const GenericMatrix>) System preconditioner matrix. nonlinear_problem (NonlinearProblem ) The nonlinear problem. iteration (std::size_t) Newton iteration number.

Parameters:
  • A
  • P
  • nonlinear_problem
  • iteration
void dolfin::NewtonSolver::update_solution(GenericVector &x, const GenericVector &dx, double relaxation_parameter, const NonlinearProblem &nonlinear_problem, std::size_t iteration)

Update solution vector by computed Newton step. Default update is given by formula:: x -= relaxation_parameter*dx Arguments x (GenericVector >) The solution vector to be updated. dx (GenericVector >) The update vector computed by Newton step. relaxation_parameter (double) Newton relaxation parameter. nonlinear_problem (NonlinearProblem ) The nonlinear problem. iteration (std::size_t) Newton iteration number.

Parameters:
  • x
  • dx
  • relaxation_parameter
  • nonlinear_problem
  • iteration
dolfin::NewtonSolver::~NewtonSolver()

Destructor.

NonlinearProblem

C++ documentation for NonlinearProblem from dolfin/nls/NonlinearProblem.h:

class dolfin::NonlinearProblem

This is a base class for nonlinear problems which can return the nonlinear function F(u) and its Jacobian J = dF(u)/du.

void dolfin::NonlinearProblem::F(GenericVector &b, const GenericVector &x) = 0

Compute F at current point x.

Parameters:
  • b
  • x
void dolfin::NonlinearProblem::J(GenericMatrix &A, const GenericVector &x) = 0

Compute J = F’ at current point x.

Parameters:
  • A
  • x
void dolfin::NonlinearProblem::J_pc(GenericMatrix &P, const GenericVector &x)

Compute J_pc used to precondition J. Not implementing this or leaving P empty results in system matrix A being used to construct preconditioner. Note that if nonempty P is not assembled on first call then a solver implementation may throw away P and not call this routine ever again.

Parameters:
  • P
  • x
dolfin::NonlinearProblem::NonlinearProblem()

Constructor.

void dolfin::NonlinearProblem::form(GenericMatrix &A, GenericMatrix &P, GenericVector &b, const GenericVector &x)

Function called by Newton solver before requesting F, J or J_pc. This can be used to compute F, J and J_pc together. Preconditioner matrix P can be left empty so that A is used instead

Parameters:
  • A
  • P
  • b
  • x
void dolfin::NonlinearProblem::form(GenericMatrix &A, GenericVector &b, const GenericVector &x)

Function called by Newton solver before requesting F or J. This can be used to compute F and J together. NOTE: This function is deprecated. Use variant with preconditioner

Parameters:
  • A
  • b
  • x
dolfin::NonlinearProblem::~NonlinearProblem()

Destructor.

OptimisationProblem

C++ documentation for OptimisationProblem from dolfin/nls/OptimisationProblem.h:

class dolfin::OptimisationProblem

This is a base class for nonlinear optimisation problems which return the real-valued objective function \(f(x)\) , its gradient \(F(x) = f'(x)``and its Hessian :math:`\) J(x) = f’‘(x)`

void dolfin::OptimisationProblem::F(GenericVector &b, const GenericVector &x) = 0

Compute the gradient \(F(x) = f'(x)\).

Parameters:
  • b
  • x
void dolfin::OptimisationProblem::J(GenericMatrix &A, const GenericVector &x) = 0

Compute the Hessian \(J(x) = f''(x)\).

Parameters:
  • A
  • x
void dolfin::OptimisationProblem::J_pc(GenericMatrix &P, const GenericVector &x)

Compute J_pc used to precondition J. Not implementing this or leaving P empty results in system matrix A being used to construct preconditioner. Note that if nonempty P is not assembled on first call then a solver implementation may throw away P and not call this routine ever again.

Parameters:
  • P
  • x
dolfin::OptimisationProblem::OptimisationProblem()

Constructor.

double dolfin::OptimisationProblem::f(const GenericVector &x) = 0

Compute the objective function \(f(x)\)

Parameters:x
void dolfin::OptimisationProblem::form(GenericMatrix &A, GenericMatrix &P, GenericVector &b, const GenericVector &x)

Function called by the solver before requesting F, J or J_pc. This can be used to compute F, J and J_pc together. Preconditioner matrix P can be left empty so that A is used instead

Parameters:
  • A
  • P
  • b
  • x
void dolfin::OptimisationProblem::form(GenericMatrix &A, GenericVector &b, const GenericVector &x)

Compute the Hessian \(J(x)=f''(x)``and the gradient :math:`\) F(x)=f’(x)` together. Called before requesting F or J. NOTE: This function is deprecated. Use variant with preconditioner

Parameters:
  • A
  • b
  • x
dolfin::OptimisationProblem::~OptimisationProblem()

Destructor.

PETScSNESSolver

C++ documentation for PETScSNESSolver from dolfin/nls/PETScSNESSolver.h:

class dolfin::PETScSNESSolver

This class implements methods for solving nonlinear systems via PETSc’s SNES interface. It includes line search and trust region techniques for globalising the convergence of the nonlinear iteration.

PetscErrorCode dolfin::PETScSNESSolver::FormFunction(SNES snes, Vec x, Vec f, void *ctx)
Parameters:
  • snes
  • x
  • f
  • ctx
PetscErrorCode dolfin::PETScSNESSolver::FormJacobian(SNES snes, Vec x, Mat A, Mat B, void *ctx)
Parameters:
  • snes
  • x
  • A
  • B
  • ctx
PetscErrorCode dolfin::PETScSNESSolver::FormObjective(SNES snes, Vec x, PetscReal *out, void *ctx)
Parameters:
  • snes
  • x
  • out
  • ctx
dolfin::PETScSNESSolver::PETScSNESSolver(MPI_Comm comm, std::string nls_type = "default")

Create SNES solver.

Parameters:
  • comm
  • nls_type
dolfin::PETScSNESSolver::PETScSNESSolver(std::string nls_type = "default")

Create SNES solver on MPI_COMM_WORLD.

Parameters:nls_type
std::string dolfin::PETScSNESSolver::get_options_prefix() const

Returns the prefix used by PETSc when searching the PETSc options database

void dolfin::PETScSNESSolver::init(NonlinearProblem &nonlinear_problem, GenericVector &x)

Set up the SNES object, but don’t do anything yet, in case the user wants to access the SNES object directly

Parameters:
  • nonlinear_problem
  • x
bool dolfin::PETScSNESSolver::is_vi() const
std::shared_ptr<const PETScVector> dolfin::PETScSNESSolver::lb
std::vector<std::pair<std::string, std::string>> dolfin::PETScSNESSolver::methods()

Return a list of available solver methods.

MPI_Comm dolfin::PETScSNESSolver::mpi_comm() const

Return the MPI communicator.

Parameters dolfin::PETScSNESSolver::parameters

Parameters .

void dolfin::PETScSNESSolver::set_bounds(GenericVector &x)
Parameters:x
void dolfin::PETScSNESSolver::set_from_options() const

Set options from the PETSc options database.

void dolfin::PETScSNESSolver::set_linear_solver_parameters()
void dolfin::PETScSNESSolver::set_options_prefix(std::string options_prefix)

Sets the prefix used by PETSc when searching the PETSc options database

Parameters:options_prefix
SNES dolfin::PETScSNESSolver::snes() const

Return PETSc SNES pointer.

std::pair<std::size_t, bool> dolfin::PETScSNESSolver::solve(NonlinearProblem &nonlinear_function, GenericVector &x)

Solve abstract nonlinear problem \(F(x) = 0\) for given \(F\) and Jacobian \(\dfrac{\partial F}{\partial x}\) . Arguments nonlinear_function (NonlinearProblem ) The nonlinear problem. x (GenericVector ) The vector. Returns std::pair<std::size_t, bool> Pair of number of Newton iterations, and whether iteration converged)

Parameters:
  • nonlinear_function
  • x
std::pair<std::size_t, bool> dolfin::PETScSNESSolver::solve(NonlinearProblem &nonlinear_problem, GenericVector &x, const GenericVector &lb, const GenericVector &ub)

Solve a nonlinear variational inequality with bound constraints Arguments nonlinear_function (NonlinearProblem ) The nonlinear problem. x (GenericVector ) The vector. lb (GenericVector ) The lower bound. ub (GenericVector ) The upper bound. Returns std::pair<std::size_t, bool> Pair of number of Newton iterations, and whether iteration converged)

Parameters:
  • nonlinear_problem
  • x
  • lb
  • ub
std::shared_ptr<const PETScVector> dolfin::PETScSNESSolver::ub
dolfin::PETScSNESSolver::~PETScSNESSolver()

Destructor.

PETScTAOSolver

C++ documentation for PETScTAOSolver from dolfin/nls/PETScTAOSolver.h:

class dolfin::PETScTAOSolver

This class implements methods for solving nonlinear optimisation problems via PETSc TAO solver. It supports unconstrained as well as bound-constrained minimisation problem

PetscErrorCode dolfin::PETScTAOSolver::FormFunctionGradient(Tao tao, Vec x, PetscReal *fobj, Vec G, void *ctx)
Parameters:
  • tao
  • x
  • fobj
  • G
  • ctx
PetscErrorCode dolfin::PETScTAOSolver::FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
Parameters:
  • tao
  • x
  • H
  • Hpre
  • ctx
dolfin::PETScTAOSolver::PETScTAOSolver(MPI_Comm comm, std::string tao_type = "default", std::string ksp_type = "default", std::string pc_type = "default")

Create TAO solver.

Parameters:
  • comm
  • tao_type
  • ksp_type
  • pc_type
dolfin::PETScTAOSolver::PETScTAOSolver(std::string tao_type = "default", std::string ksp_type = "default", std::string pc_type = "default")

Create TAO solver on MPI_COMM_WORLD.

Parameters:
  • tao_type
  • ksp_type
  • pc_type
PetscErrorCode dolfin::PETScTAOSolver::TaoConvergenceTest(Tao tao, void *ctx)
Parameters:
  • tao
  • ctx
void dolfin::PETScTAOSolver::init(OptimisationProblem &optimisation_problem, PETScVector &x)

Initialise the TAO solver for an unconstrained minimisation problem, in case the user wants to access the TAO object directly

Parameters:
  • optimisation_problem
  • x
void dolfin::PETScTAOSolver::init(OptimisationProblem &optimisation_problem, PETScVector &x, const PETScVector &lb, const PETScVector &ub)

Initialise the TAO solver for a bound-constrained minimisation problem, in case the user wants to access the TAO object directly

Parameters:
  • optimisation_problem
  • x
  • lb
  • ub
std::vector<std::pair<std::string, std::string>> dolfin::PETScTAOSolver::methods()

Return a list of available solver methods.

MPI_Comm dolfin::PETScTAOSolver::mpi_comm() const

Return the MPI communicator.

Parameters dolfin::PETScTAOSolver::parameters

Parameters for the PETSc TAO solver.

void dolfin::PETScTAOSolver::set_ksp_options()
void dolfin::PETScTAOSolver::set_tao(std::string tao_type)
Parameters:tao_type
void dolfin::PETScTAOSolver::set_tao_options()
std::pair<std::size_t, bool> dolfin::PETScTAOSolver::solve(OptimisationProblem &optimisation_problem, GenericVector &x)

Solve a nonlinear unconstrained minimisation problem Arguments optimisation_problem (OptimisationProblem ) The nonlinear optimisation problem. x (GenericVector ) The solution vector (initial guess). Returns (its, converged) (std::pair<std::size_t, bool>) Pair of number of iterations, and whether iteration converged

Parameters:
  • optimisation_problem
  • x
std::pair<std::size_t, bool> dolfin::PETScTAOSolver::solve(OptimisationProblem &optimisation_problem, GenericVector &x, const GenericVector &lb, const GenericVector &ub)

Solve a nonlinear bound-constrained optimisation problem Arguments optimisation_problem (OptimisationProblem ) The nonlinear optimisation problem. x (GenericVector ) The solution vector (initial guess). lb (GenericVector ) The lower bound. ub (GenericVector ) The upper bound. Returns (its, converged) (std::pair<std::size_t, bool>) Pair of number of iterations, and whether iteration converged

Parameters:
  • optimisation_problem
  • x
  • lb
  • ub
std::pair<std::size_t, bool> dolfin::PETScTAOSolver::solve(OptimisationProblem &optimisation_problem, PETScVector &x, const PETScVector &lb, const PETScVector &ub)

Solve a nonlinear bound-constrained minimisation problem Arguments optimisation_problem (OptimisationProblem ) The nonlinear optimisation problem. x (PETScVector ) The solution vector (initial guess). lb (PETScVector ) The lower bound. ub (PETScVector ) The upper bound. Returns (its, converged) (std::pair<std::size_t, bool>) Pair of number of iterations, and whether iteration converged

Parameters:
  • optimisation_problem
  • x
  • lb
  • ub
Tao dolfin::PETScTAOSolver::tao() const

Return the TAO pointer.

void dolfin::PETScTAOSolver::update_parameters(std::string tao_type, std::string ksp_type, std::string pc_type)
Parameters:
  • tao_type
  • ksp_type
  • pc_type
dolfin::PETScTAOSolver::~PETScTAOSolver()

Destructor.

TAOLinearBoundSolver

C++ documentation for TAOLinearBoundSolver from dolfin/nls/TAOLinearBoundSolver.h:

class dolfin::TAOLinearBoundSolver

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 \(x_l\leq x\leq x_u\) such that \((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.

# 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)
dolfin::TAOLinearBoundSolver::TAOLinearBoundSolver(MPI_Comm comm)

Create TAO bound constrained solver.

Parameters:comm
dolfin::TAOLinearBoundSolver::TAOLinearBoundSolver(const std::string method = "default", const std::string ksp_type = "default", const std::string pc_type = "default")

Create TAO bound constrained solver.

Parameters:
  • method
  • ksp_type
  • pc_type
std::shared_ptr<const PETScMatrix> dolfin::TAOLinearBoundSolver::get_matrix() const

Return Matrix shared pointer.

std::shared_ptr<const PETScVector> dolfin::TAOLinearBoundSolver::get_vector() const

Return load vector shared pointer.

void dolfin::TAOLinearBoundSolver::init(const std::string &method)
Parameters:method
std::map<std::string, std::string> dolfin::TAOLinearBoundSolver::krylov_solvers()

Return a list of available krylov solvers.

std::map<std::string, std::string> dolfin::TAOLinearBoundSolver::methods()

Return a list of available Tao solver methods.

std::map<std::string, std::string> dolfin::TAOLinearBoundSolver::preconditioners()

Return a list of available preconditioners.

void dolfin::TAOLinearBoundSolver::read_parameters()
void dolfin::TAOLinearBoundSolver::set_ksp(const std::string ksp_type = "default")

Set PETSC Krylov Solver (ksp) used by TAO.

Parameters:ksp_type
void dolfin::TAOLinearBoundSolver::set_ksp_options()
void dolfin::TAOLinearBoundSolver::set_operators(std::shared_ptr<const GenericMatrix> A, std::shared_ptr<const GenericVector> b)
Parameters:
  • A
  • b
void dolfin::TAOLinearBoundSolver::set_operators(std::shared_ptr<const PETScMatrix> A, std::shared_ptr<const PETScVector> b)
Parameters:
  • A
  • b
void dolfin::TAOLinearBoundSolver::set_solver(const std::string&)

Set the TAO solver type.

Parameters:method
std::size_t dolfin::TAOLinearBoundSolver::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

Parameters:
  • A
  • x
  • b
  • xl
  • xu
std::size_t dolfin::TAOLinearBoundSolver::solve(const PETScMatrix &A, PETScVector &x, const PETScVector &b, const PETScVector &xl, const PETScVector &xu)

Solve the linear variational inequality defined by A and b with xl =< x <= xu

Parameters:
  • A
  • x
  • b
  • xl
  • xu
Tao dolfin::TAOLinearBoundSolver::tao() const

Return TAO solver pointer.

dolfin::TAOLinearBoundSolver::~TAOLinearBoundSolver()

Destructor.