Class implementing LBFGS optimzation method.
More...
#include <LBFGSNonLinearSolver.h>
|
| LBFGSNonLinearSolver (const bool usePreconditioner, const double maxUpdate, const unsigned int maxNumberIterations, const int maxNumPastSteps, const unsigned int debugLevel, const MPI_Comm &mpi_comm_parent, const bool isCurvatureOnlyLineSearchStoppingCondition=false) |
| Constructor.
|
|
| ~LBFGSNonLinearSolver () |
| Destructor.
|
|
nonLinearSolver::ReturnValueType | solve (nonlinearSolverProblem &problem, const std::string checkpointFileName="", const bool restart=false) |
| Solve non-linear problem using LBFGS method.
|
|
void | save (const std::string &checkpointFileName) |
| Create checkpoint file for current state of the LBFGS solver.
|
|
virtual | ~nonLinearSolver ()=0 |
| Destructor.
|
|
Class implementing LBFGS optimzation method.
- Author
- Nikhil Kodali
◆ LBFGSNonLinearSolver()
dftfe::LBFGSNonLinearSolver::LBFGSNonLinearSolver |
( |
const bool | usePreconditioner, |
|
|
const double | maxUpdate, |
|
|
const unsigned int | maxNumberIterations, |
|
|
const int | maxNumPastSteps, |
|
|
const unsigned int | debugLevel, |
|
|
const MPI_Comm & | mpi_comm_parent, |
|
|
const bool | isCurvatureOnlyLineSearchStoppingCondition = false ) |
Constructor.
- Parameters
-
usePreconditioner | Boolean parameter specifying whether or not to use the preconditioner. |
maxUpdate | Maximum allowed step length in any direction. |
maxNumberIterations | Maximum number of iterations. |
maxNumPastSteps | Number of previous steps stored by the LBFGS solver. |
debugLevel | Debug output level: 0 - no debug output 1 - limited debug output 2 - all debug output. |
mpi_comm_parent | The mpi communicator used. |
◆ ~LBFGSNonLinearSolver()
dftfe::LBFGSNonLinearSolver::~LBFGSNonLinearSolver |
( |
| ) |
|
◆ checkWolfe()
void dftfe::LBFGSNonLinearSolver::checkWolfe |
( |
| ) |
|
|
private |
Test if the step satisfies strong Wolfe conditions.
◆ computeHx()
void dftfe::LBFGSNonLinearSolver::computeHx |
( |
std::vector< double > & | Hx | ) |
|
|
private |
Compute Hessian inverse times vector.
◆ computeStep()
void dftfe::LBFGSNonLinearSolver::computeStep |
( |
| ) |
|
|
private |
◆ computeStepScale()
Compute scaling factor for the step.
◆ computeUpdateStep()
void dftfe::LBFGSNonLinearSolver::computeUpdateStep |
( |
| ) |
|
|
private |
◆ initializePreconditioner()
Initialize preconditioner.
◆ load()
void dftfe::LBFGSNonLinearSolver::load |
( |
const std::string & | checkpointFileName | ) |
|
|
private |
Load LBFGS solver state from checkpoint file.
◆ save()
void dftfe::LBFGSNonLinearSolver::save |
( |
const std::string & | checkpointFileName | ) |
|
|
virtual |
◆ scalePreconditioner()
◆ solve()
Solve non-linear problem using LBFGS 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.
◆ updateHistory()
void dftfe::LBFGSNonLinearSolver::updateHistory |
( |
| ) |
|
|
private |
Update the stored history, damped LBFGS.
◆ updateSolution()
bool dftfe::LBFGSNonLinearSolver::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_alpha
double dftfe::LBFGSNonLinearSolver::d_alpha |
|
private |
storage for backtracking line search parameter.
◆ d_deltaGq
std::deque<std::vector<double> > dftfe::LBFGSNonLinearSolver::d_deltaGq |
|
private |
◆ d_deltaX
std::vector<double> dftfe::LBFGSNonLinearSolver::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::LBFGSNonLinearSolver::d_deltaXNew |
|
private |
◆ d_deltaXq
std::deque<std::vector<double> > dftfe::LBFGSNonLinearSolver::d_deltaXq |
|
private |
◆ d_gradient
std::vector<double> dftfe::LBFGSNonLinearSolver::d_gradient |
|
private |
storage for the value and gradient of the nonlinear problem in the current bfgs step
◆ d_gradientNew
std::vector<double> dftfe::LBFGSNonLinearSolver::d_gradientNew |
|
private |
storage for the value and gradient of the nonlinear problem evaluated at the end of the current bfgs step
◆ d_isCurvatureOnlyLineSearchStoppingCondition
bool dftfe::LBFGSNonLinearSolver::d_isCurvatureOnlyLineSearchStoppingCondition |
|
private |
◆ d_iter
unsigned int dftfe::LBFGSNonLinearSolver::d_iter |
|
private |
storage for current bfgs iteration count.
◆ d_maxNumPastSteps
const int dftfe::LBFGSNonLinearSolver::d_maxNumPastSteps |
|
private |
storage for the maximum number of past steps to be stored.
◆ d_maxStepLength
double dftfe::LBFGSNonLinearSolver::d_maxStepLength |
|
private |
storage for the maximum allowed step length.
◆ d_noHistory
bool dftfe::LBFGSNonLinearSolver::d_noHistory |
|
private |
◆ d_normDeltaXnew
double dftfe::LBFGSNonLinearSolver::d_normDeltaXnew |
|
private |
storage for inf norm of the update step.
◆ d_numberUnknowns
unsigned int dftfe::LBFGSNonLinearSolver::d_numberUnknowns |
|
private |
storage for number of unknowns to be solved for in the nonlinear problem.
◆ d_numPastSteps
int dftfe::LBFGSNonLinearSolver::d_numPastSteps |
|
private |
storage for the number of past steps currently stored.
◆ d_preconditioner
std::vector<double> dftfe::LBFGSNonLinearSolver::d_preconditioner |
|
private |
storage for the preconditioner.
◆ d_rhoq
std::deque<double> dftfe::LBFGSNonLinearSolver::d_rhoq |
|
private |
◆ d_stepAccepted
bool dftfe::LBFGSNonLinearSolver::d_stepAccepted |
|
private |
boolean parameter for step accepteance and Wolfe conditions.
◆ d_updateVector
std::vector<double> dftfe::LBFGSNonLinearSolver::d_updateVector |
|
private |
◆ d_usePreconditioner
const bool dftfe::LBFGSNonLinearSolver::d_usePreconditioner |
|
private |
flag for using the preconditioner
◆ d_useSingleAtomSolutionsInitialGuess
bool dftfe::LBFGSNonLinearSolver::d_useSingleAtomSolutionsInitialGuess |
|
private |
◆ d_value
std::vector<double> dftfe::LBFGSNonLinearSolver::d_value |
|
private |
◆ d_valueNew
std::vector<double> dftfe::LBFGSNonLinearSolver::d_valueNew |
|
private |
◆ d_wolfeCurvature
bool dftfe::LBFGSNonLinearSolver::d_wolfeCurvature |
|
private |
◆ d_wolfeSatisfied
bool dftfe::LBFGSNonLinearSolver::d_wolfeSatisfied |
|
private |
◆ d_wolfeSufficientDec
bool dftfe::LBFGSNonLinearSolver::d_wolfeSufficientDec |
|
private |
◆ mpi_communicator
MPI_Comm dftfe::LBFGSNonLinearSolver::mpi_communicator |
|
private |
◆ pcout
dealii::ConditionalOStream dftfe::LBFGSNonLinearSolver::pcout |
|
private |
The documentation for this class was generated from the following file: