DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
LBFGSNonLinearSolver.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17
18#ifndef LBFGSNonLinearSolver_h
19#define LBFGSNonLinearSolver_h
20
21
22#include "nonLinearSolver.h"
24namespace dftfe
25{
26 /**
27 * @brief Class implementing LBFGS optimzation method.
28 *
29 * @author Nikhil Kodali
30 */
32 {
33 public:
34 /**
35 * @brief Constructor.
36 *
37 * @param usePreconditioner Boolean parameter specifying whether or not to use the preconditioner.
38 * @param maxUpdate Maximum allowed step length in any direction.
39 * @param maxNumberIterations Maximum number of iterations.
40 * @param maxNumPastSteps Number of previous steps stored by the LBFGS solver.
41 * @param debugLevel Debug output level:
42 * 0 - no debug output
43 * 1 - limited debug output
44 * 2 - all debug output.
45 * @param mpi_comm_parent The mpi communicator used.
46 */
48 const bool usePreconditioner,
49 const double maxUpdate,
50 const unsigned int maxNumberIterations,
51 const int maxNumPastSteps,
52 const unsigned int debugLevel,
53 const MPI_Comm & mpi_comm_parent,
54 const bool isCurvatureOnlyLineSearchStoppingCondition = false);
55
56 /**
57 * @brief Destructor.
58 */
60
61 /**
62 * @brief Solve non-linear problem using LBFGS method.
63 *
64 * @param problem[in] nonlinearSolverProblem object.
65 * @param checkpointFileName[in] if string is non-empty, creates checkpoint file
66 * named checkpointFileName for every nonlinear iteration. If restart is set
67 * to true, checkpointFileName must match the name of the checkpoint file.
68 * Empty string will throw an error.
69 * @param restart[in] boolean specifying whether this is a restart solve using the checkpoint file
70 * specified by checkpointFileName.
71 * @return Return value indicating success or failure.
72 */
75 const std::string checkpointFileName = "",
76 const bool restart = false);
77
78
79
80 /**
81 * @brief Create checkpoint file for current state of the LBFGS solver.
82 *
83 */
84 void
85 save(const std::string &checkpointFileName);
86
87
88 private:
89 /**
90 * @brief Initialize preconditioner.
91 */
92 void
94
95 /**
96 * @brief Scale preconditioner.
97 */
98 void
100
101 /**
102 * @brief Compute Hessian inverse times vector.
103 */
104 void
105 computeHx(std::vector<double> &Hx);
106
107 /**
108 * @brief Compute LBFGS step.
109 */
110 void
112
113 /**
114 * @brief Compute Update Vector.
115 */
116 void
118
119 /**
120 * @brief Update the stored history, damped LBFGS.
121 */
122 void
124
125 /**
126 * @brief Test if the step satisfies strong Wolfe conditions.
127 */
128 void
130
131 /**
132 * @brief Compute scaling factor for the step.
133 */
134 void
136
137 /**
138 * @brief Update solution x -> x + step.
139 *
140 * @param step update step vector.
141 * @param problem nonlinearSolverProblem object.
142 *
143 * @return bool true if valid update and false if increment bound exceeded
144 *
145 */
146 bool
147 updateSolution(const std::vector<double> &step,
148 nonlinearSolverProblem & problem);
149
150 /**
151 * @brief Load LBFGS solver state from checkpoint file.
152 *
153 */
154 void
155 load(const std::string &checkpointFileName);
156
157 /// storage for the value and gradient of the nonlinear problem in the
158 /// current bfgs step
159 std::vector<double> d_gradient, d_value;
160
161 /// storage for the value and gradient of the nonlinear problem evaluated at
162 /// the end of the current bfgs step
163 std::vector<double> d_gradientNew, d_valueNew;
164
165 /// storage for the step taken in last BFGS step, step computed in the
166 /// corrent BFGS step and the update vector computed in the current bfgs
167 /// step.
168 std::vector<double> d_deltaX, d_deltaXNew, d_updateVector;
169
170 /// storage for the preconditioner.
171 std::vector<double> d_preconditioner;
172
173 /// storage for number of unknowns to be solved for in the nonlinear
174 /// problem.
175 unsigned int d_numberUnknowns;
176
177 /// storage for current bfgs iteration count.
178 unsigned int d_iter;
179
180 /// storage for history.
181 std::deque<std::vector<double>> d_deltaGq, d_deltaXq;
182 std::deque<double> d_rhoq;
183
184 /// storage for the maximum number of past steps to be stored.
186
187 /// storage for the number of past steps currently stored.
189
190 /// storage for inf norm of the update step.
192
193 /// storage for the maximum allowed step length.
195
196 /// storage for backtracking line search parameter.
197 double d_alpha;
198
199 /// boolean parameter for step accepteance and Wolfe conditions.
202
203 /// flag for using the preconditioner
205
207
208 //
210
211
212
213 // parallel objects
215 dealii::ConditionalOStream pcout;
216 };
217
218} // namespace dftfe
219#endif // BFGSNonLinearSolver_h
void computeStepScale(nonlinearSolverProblem &problem)
Compute scaling factor for the step.
bool d_wolfeSufficientDec
Definition LBFGSNonLinearSolver.h:200
std::vector< double > d_preconditioner
storage for the preconditioner.
Definition LBFGSNonLinearSolver.h:171
void updateHistory()
Update the stored history, damped LBFGS.
void save(const std::string &checkpointFileName)
Create checkpoint file for current state of the LBFGS solver.
void initializePreconditioner(nonlinearSolverProblem &problem)
Initialize preconditioner.
void load(const std::string &checkpointFileName)
Load LBFGS solver state from checkpoint file.
dealii::ConditionalOStream pcout
Definition LBFGSNonLinearSolver.h:215
void computeUpdateStep()
Compute Update Vector.
std::vector< double > d_gradientNew
Definition LBFGSNonLinearSolver.h:163
bool d_useSingleAtomSolutionsInitialGuess
Definition LBFGSNonLinearSolver.h:206
void computeStep()
Compute LBFGS step.
std::vector< double > d_updateVector
Definition LBFGSNonLinearSolver.h:168
double d_normDeltaXnew
storage for inf norm of the update step.
Definition LBFGSNonLinearSolver.h:191
std::vector< double > d_value
Definition LBFGSNonLinearSolver.h:159
double d_alpha
storage for backtracking line search parameter.
Definition LBFGSNonLinearSolver.h:197
std::vector< double > d_gradient
Definition LBFGSNonLinearSolver.h:159
bool d_wolfeSatisfied
Definition LBFGSNonLinearSolver.h:201
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.
const bool d_usePreconditioner
flag for using the preconditioner
Definition LBFGSNonLinearSolver.h:204
int d_numPastSteps
storage for the number of past steps currently stored.
Definition LBFGSNonLinearSolver.h:188
void scalePreconditioner(nonlinearSolverProblem &problem)
Scale preconditioner.
unsigned int d_numberUnknowns
Definition LBFGSNonLinearSolver.h:175
const int d_maxNumPastSteps
storage for the maximum number of past steps to be stored.
Definition LBFGSNonLinearSolver.h:185
std::deque< std::vector< double > > d_deltaGq
storage for history.
Definition LBFGSNonLinearSolver.h:181
unsigned int d_iter
storage for current bfgs iteration count.
Definition LBFGSNonLinearSolver.h:178
std::vector< double > d_deltaX
Definition LBFGSNonLinearSolver.h:168
void computeHx(std::vector< double > &Hx)
Compute Hessian inverse times vector.
std::vector< double > d_deltaXNew
Definition LBFGSNonLinearSolver.h:168
bool d_stepAccepted
boolean parameter for step accepteance and Wolfe conditions.
Definition LBFGSNonLinearSolver.h:200
bool d_isCurvatureOnlyLineSearchStoppingCondition
Definition LBFGSNonLinearSolver.h:209
~LBFGSNonLinearSolver()
Destructor.
bool updateSolution(const std::vector< double > &step, nonlinearSolverProblem &problem)
Update solution x -> x + step.
std::vector< double > d_valueNew
Definition LBFGSNonLinearSolver.h:163
void checkWolfe()
Test if the step satisfies strong Wolfe conditions.
bool d_noHistory
Definition LBFGSNonLinearSolver.h:206
std::deque< std::vector< double > > d_deltaXq
Definition LBFGSNonLinearSolver.h:181
bool d_wolfeCurvature
Definition LBFGSNonLinearSolver.h:200
std::deque< double > d_rhoq
Definition LBFGSNonLinearSolver.h:182
nonLinearSolver::ReturnValueType solve(nonlinearSolverProblem &problem, const std::string checkpointFileName="", const bool restart=false)
Solve non-linear problem using LBFGS method.
double d_maxStepLength
storage for the maximum allowed step length.
Definition LBFGSNonLinearSolver.h:194
MPI_Comm mpi_communicator
Definition LBFGSNonLinearSolver.h:214
nonLinearSolver(const unsigned int debugLevel, const unsigned int maxNumberIterations, const double tolerance=0)
Constructor.
ReturnValueType
Definition nonLinearSolver.h:45
Abstract class for solver functions.
Definition nonlinearSolverProblem.h:30
Definition pseudoPotentialToDftfeConverter.cc:34