Class implementing a modified BFGS optimization scheme.
More...
#include <BFGSNonLinearSolver.h>
|
| BFGSNonLinearSolver (const bool usePreconditioner, const bool useRFOStep, const unsigned int maxNumberIterations, const unsigned int debugLevel, const MPI_Comm &mpi_comm_parent, const double trustRadius_maximum=0.5, const double trustRadius_initial=0.02, const double trustRadius_minimum=1.0e-4, const bool isCurvatureOnlyLineSearchStoppingCondition=false) |
| Constructor.
|
|
| ~BFGSNonLinearSolver () |
| Destructor.
|
|
nonLinearSolver::ReturnValueType | solve (nonlinearSolverProblem &problem, const std::string checkpointFileName="", const bool restart=false) |
| Solve non-linear problem using a modified BFGS method.
|
|
void | save (const std::string &checkpointFileName) |
| Create checkpoint file for current state of the BFGS solver.
|
|
virtual | ~nonLinearSolver ()=0 |
| Destructor.
|
|
|
void | initializeHessian (nonlinearSolverProblem &problem) |
| initialize hessian, either preconditioner or identity matrix.
|
|
void | updateHessian () |
| Update Hessian according to damped BFGS rule: Procedure 18.2 of Nocedal and Wright.
|
|
void | scaleHessian () |
| Scale hessian according to eqn 6.20 of Nocedal and Wright.
|
|
void | checkWolfe () |
| Check if the step satifies the Strong Wolfe conditons.
|
|
void | computeRFOStep () |
| Compute step using the Rational Function Method.
|
|
void | computeNewtonStep () |
| Compute the Quasi-Newton Step.
|
|
void | computeStep () |
| Compute the final update step using the trust radius and whether or not the previous step was accepted.
|
|
void | computeTrustRadius (nonlinearSolverProblem &problem) |
| Estimate the trust radius for the next step based on the previous step and check for trust radius max/min conditons and reset BFGS if needed.
|
|
bool | updateSolution (const std::vector< double > &step, nonlinearSolverProblem &problem) |
| Update solution x -> x + step.
|
|
void | load (const std::string &checkpointFileName) |
| Load BFGS solver state from checkpoint file.
|
|
Class implementing a modified BFGS optimization scheme.
- Author
- Nikhil Kodali
◆ BFGSNonLinearSolver()
dftfe::BFGSNonLinearSolver::BFGSNonLinearSolver |
( |
const bool | usePreconditioner, |
|
|
const bool | useRFOStep, |
|
|
const unsigned int | maxNumberIterations, |
|
|
const unsigned int | debugLevel, |
|
|
const MPI_Comm & | mpi_comm_parent, |
|
|
const double | trustRadius_maximum = 0.5, |
|
|
const double | trustRadius_initial = 0.02, |
|
|
const double | trustRadius_minimum = 1.0e-4, |
|
|
const bool | isCurvatureOnlyLineSearchStoppingCondition = false ) |
Constructor.
- Parameters
-
usePreconditioner | Boolean parameter specifying whether or not to use the preconditioner. |
useRFOStep | Boolean parameter specifying whether or not the RFO step is used. |
maxNumberIterations | Maximum number of iterations. |
debugLevel | Debug output level: 0 - no debug output 1 - limited debug output 2 - all debug output. |
mpi_comm_parent | The mpi communicator used. |
trustRadius_maximum | Maximum trust region radius. |
trustRadius_initial | Initial trust region radius. |
trustRadius_minimum | mimimum trust region radius (will reset BFGS). |
◆ ~BFGSNonLinearSolver()
dftfe::BFGSNonLinearSolver::~BFGSNonLinearSolver |
( |
| ) |
|
◆ checkWolfe()
void dftfe::BFGSNonLinearSolver::checkWolfe |
( |
| ) |
|
|
private |
Check if the step satifies the Strong Wolfe conditons.
◆ computeNewtonStep()
void dftfe::BFGSNonLinearSolver::computeNewtonStep |
( |
| ) |
|
|
private |
Compute the Quasi-Newton Step.
◆ computeRFOStep()
void dftfe::BFGSNonLinearSolver::computeRFOStep |
( |
| ) |
|
|
private |
Compute step using the Rational Function Method.
◆ computeStep()
void dftfe::BFGSNonLinearSolver::computeStep |
( |
| ) |
|
|
private |
Compute the final update step using the trust radius and whether or not the previous step was accepted.
◆ computeTrustRadius()
Estimate the trust radius for the next step based on the previous step and check for trust radius max/min conditons and reset BFGS if needed.
◆ initializeHessian()
initialize hessian, either preconditioner or identity matrix.
◆ load()
void dftfe::BFGSNonLinearSolver::load |
( |
const std::string & | checkpointFileName | ) |
|
|
private |
Load BFGS solver state from checkpoint file.
◆ save()
void dftfe::BFGSNonLinearSolver::save |
( |
const std::string & | checkpointFileName | ) |
|
|
virtual |
◆ scaleHessian()
void dftfe::BFGSNonLinearSolver::scaleHessian |
( |
| ) |
|
|
private |
Scale hessian according to eqn 6.20 of Nocedal and Wright.
◆ solve()
Solve non-linear problem using a modified BFGS method.
- Parameters
-
problem[in] | nonlinearSolverProblem object. |
checkpointFileName[in] | if string is non-empty, creates checkpoint file named checkpointFileName for every nonlinear iteration. If restart is set to true, checkpointFileName must match the name of the checkpoint file. Empty string will throw an error. |
restart[in] | boolean specifying whether this is a restart solve using the checkpoint file specified by checkpointFileName. |
- Returns
- Return value indicating success or failure.
Implements dftfe::nonLinearSolver.
◆ updateHessian()
void dftfe::BFGSNonLinearSolver::updateHessian |
( |
| ) |
|
|
private |
Update Hessian according to damped BFGS rule: Procedure 18.2 of Nocedal and Wright.
◆ updateSolution()
bool dftfe::BFGSNonLinearSolver::updateSolution |
( |
const std::vector< double > & | step, |
|
|
nonlinearSolverProblem & | problem ) |
|
private |
Update solution x -> x + step.
- Parameters
-
step | update step vector. |
problem | nonlinearSolverProblem object. |
- Returns
- bool true if valid update and false if increment bound exceeded
◆ d_deltaX
std::vector<double> dftfe::BFGSNonLinearSolver::d_deltaX |
|
private |
storage for the step taken in last BFGS step, step computed in the corrent BFGS step and the update vector computed in the current bfgs step.
◆ d_deltaXNew
std::vector<double> dftfe::BFGSNonLinearSolver::d_deltaXNew |
|
private |
◆ d_gradient
std::vector<double> dftfe::BFGSNonLinearSolver::d_gradient |
|
private |
storage for the value and gradient of the nonlinear problem in the current bfgs step.
◆ d_gradientNew
std::vector<double> dftfe::BFGSNonLinearSolver::d_gradientNew |
|
private |
storage for the value and gradient of the nonlinear problem evaluated at the end of the current bfgs step.
◆ d_hessian
std::vector<double> dftfe::BFGSNonLinearSolver::d_hessian |
|
private |
storage for the hessian in current bfgs step.
◆ d_hessianScaled
bool dftfe::BFGSNonLinearSolver::d_hessianScaled |
|
private |
flag to check if hessian is scaled.
◆ d_isCurvatureOnlyLineSearchStoppingCondition
bool dftfe::BFGSNonLinearSolver::d_isCurvatureOnlyLineSearchStoppingCondition |
|
private |
◆ d_isReset
int dftfe::BFGSNonLinearSolver::d_isReset |
|
private |
Flag to store the reset state, 0 if step is accepted, 1 if reset occured and no steps are accepted, 2 if two resets occur without step being accepted (failure of BFGS).
◆ d_iter
unsigned int dftfe::BFGSNonLinearSolver::d_iter |
|
private |
storage for current bfgs iteration count
◆ d_normDeltaXnew
double dftfe::BFGSNonLinearSolver::d_normDeltaXnew |
|
private |
storage for inf norm of the update step.
◆ d_numberUnknowns
unsigned int dftfe::BFGSNonLinearSolver::d_numberUnknowns |
|
private |
storage for number of unknowns to be solved for in the nonlinear problem.
◆ d_Srfo
std::vector<double> dftfe::BFGSNonLinearSolver::d_Srfo |
|
private |
storage for the S matrix in RFO framework, initialized to starting. Hessian Guess
◆ d_stepAccepted
bool dftfe::BFGSNonLinearSolver::d_stepAccepted |
|
private |
boolean parameter for step accepteance and Wolfe conditions.
◆ d_trustRadius
double dftfe::BFGSNonLinearSolver::d_trustRadius |
|
private |
◆ d_trustRadiusInitial
double dftfe::BFGSNonLinearSolver::d_trustRadiusInitial |
|
private |
storage for trust region parameters.
◆ d_trustRadiusMax
double dftfe::BFGSNonLinearSolver::d_trustRadiusMax |
|
private |
◆ d_trustRadiusMin
double dftfe::BFGSNonLinearSolver::d_trustRadiusMin |
|
private |
◆ d_updateVector
std::vector<double> dftfe::BFGSNonLinearSolver::d_updateVector |
|
private |
◆ d_usePreconditioner
const bool dftfe::BFGSNonLinearSolver::d_usePreconditioner |
|
private |
◆ d_useRFOStep
const bool dftfe::BFGSNonLinearSolver::d_useRFOStep |
|
private |
◆ d_useSingleAtomSolutionsInitialGuess
bool dftfe::BFGSNonLinearSolver::d_useSingleAtomSolutionsInitialGuess |
|
private |
◆ d_value
std::vector<double> dftfe::BFGSNonLinearSolver::d_value |
|
private |
◆ d_valueNew
std::vector<double> dftfe::BFGSNonLinearSolver::d_valueNew |
|
private |
◆ d_wolfeCurvature
bool dftfe::BFGSNonLinearSolver::d_wolfeCurvature |
|
private |
◆ d_wolfeSatisfied
bool dftfe::BFGSNonLinearSolver::d_wolfeSatisfied |
|
private |
◆ d_wolfeSufficientDec
bool dftfe::BFGSNonLinearSolver::d_wolfeSufficientDec |
|
private |
◆ mpi_communicator
MPI_Comm dftfe::BFGSNonLinearSolver::mpi_communicator |
|
private |
◆ pcout
dealii::ConditionalOStream dftfe::BFGSNonLinearSolver::pcout |
|
private |
The documentation for this class was generated from the following file: