DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
BFGSNonLinearSolver.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 BFGSNonLinearSolver_h
19#define BFGSNonLinearSolver_h
20
21
22#include "nonLinearSolver.h"
24namespace dftfe
25{
26 /**
27 * @brief Class implementing a modified BFGS optimization scheme
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 useRFOStep Boolean parameter specifying whether or not the RFO step is used.
39 * @param maxNumberIterations Maximum number of iterations.
40 * @param debugLevel Debug output level:
41 * 0 - no debug output
42 * 1 - limited debug output
43 * 2 - all debug output.
44 * @param mpi_comm_parent The mpi communicator used.
45 * @param trustRadius_maximum Maximum trust region radius.
46 * @param trustRadius_initial Initial trust region radius.
47 * @param trustRadius_minimum mimimum trust region radius (will reset BFGS).
48 */
50 const bool usePreconditioner,
51 const bool useRFOStep,
52 const unsigned int maxNumberIterations,
53 const unsigned int debugLevel,
54 const MPI_Comm & mpi_comm_parent,
55 const double trustRadius_maximum = 0.5,
56 const double trustRadius_initial = 0.02,
57 const double trustRadius_minimum = 1.0e-4,
58 const bool isCurvatureOnlyLineSearchStoppingCondition = false);
59
60 /**
61 * @brief Destructor.
62 */
64
65 /**
66 * @brief Solve non-linear problem using a modified BFGS method.
67 *
68 * @param problem[in] nonlinearSolverProblem object.
69 * @param checkpointFileName[in] if string is non-empty, creates checkpoint file
70 * named checkpointFileName for every nonlinear iteration. If restart is set
71 * to true, checkpointFileName must match the name of the checkpoint file.
72 * Empty string will throw an error.
73 * @param restart[in] boolean specifying whether this is a restart solve using the checkpoint file
74 * specified by checkpointFileName.
75 * @return Return value indicating success or failure.
76 */
79 const std::string checkpointFileName = "",
80 const bool restart = false);
81
82
83 /**
84 * @brief Create checkpoint file for current state of the BFGS solver.
85 *
86 */
87 void
88 save(const std::string &checkpointFileName);
89
90 private:
91 /**
92 * @brief initialize hessian, either preconditioner or identity matrix.
93 */
94 void
96
97 /**
98 * @brief Update Hessian according to damped BFGS rule: Procedure 18.2 of Nocedal and Wright.
99 */
100 void
102
103 /**
104 * @brief Scale hessian according to eqn 6.20 of Nocedal and Wright.
105 */
106 void
108
109 /**
110 * @brief Check if the step satifies the Strong Wolfe conditons.
111 */
112 void
114
115 /**
116 * @brief Compute step using the Rational Function Method.
117 */
118 void
120
121 /**
122 * @brief Compute the Quasi-Newton Step.
123 */
124 void
126
127 /**
128 * @brief Compute the final update step using the trust radius and whether or not the previous step was accepted.
129 */
130 void
132
133 /**
134 * @brief 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.
135 */
136 void
138
139
140 /**
141 * @brief Update solution x -> x + step.
142 *
143 * @param step update step vector.
144 * @param problem nonlinearSolverProblem object.
145 *
146 * @return bool true if valid update and false if increment bound exceeded
147 *
148 */
149 bool
150 updateSolution(const std::vector<double> &step,
151 nonlinearSolverProblem & problem);
152
153
154 /**
155 * @brief Load BFGS solver state from checkpoint file.
156 *
157 */
158 void
159 load(const std::string &checkpointFileName);
160
161 /// storage for the value and gradient of the nonlinear problem in the
162 /// current bfgs step.
163 std::vector<double> d_gradient, d_value;
164
165 /// storage for the value and gradient of the nonlinear problem evaluated at
166 /// the end of the current bfgs step.
167 std::vector<double> d_gradientNew, d_valueNew;
168
169 /// storage for the step taken in last BFGS step, step computed in the
170 /// corrent BFGS step and the update vector computed in the current bfgs
171 /// step.
172 std::vector<double> d_deltaX, d_deltaXNew, d_updateVector;
173
174 /// storage for number of unknowns to be solved for in the nonlinear
175 /// problem.
176 unsigned int d_numberUnknowns;
177
178 /// storage for current bfgs iteration count
179 unsigned int d_iter;
180
181 /// storage for the S matrix in RFO framework, initialized to starting.
182 /// Hessian Guess
183 std::vector<double> d_Srfo;
184
185 /// storage for the hessian in current bfgs step.
186 std::vector<double> d_hessian;
187
188 /// storage for inf norm of the update step.
190
191 /// storage for trust region parameters.
194
195 /// boolean parameter for step accepteance and Wolfe conditions.
198
199 ///
200 /// flag to check if hessian is scaled.
201 ///
203
204 //
206
207 /// Flag to store the reset state, 0 if step is accepted, 1 if reset occured
208 /// and no steps are accepted, 2 if two resets occur without step being
209 /// accepted (failure of BFGS).
211
214
215 // parallel objects
217 dealii::ConditionalOStream pcout;
218 };
219
220} // namespace dftfe
221#endif // BFGSNonLinearSolver_h
dealii::ConditionalOStream pcout
Definition BFGSNonLinearSolver.h:217
~BFGSNonLinearSolver()
Destructor.
void checkWolfe()
Check if the step satifies the Strong Wolfe conditons.
const bool d_useRFOStep
Definition BFGSNonLinearSolver.h:213
bool d_wolfeCurvature
Definition BFGSNonLinearSolver.h:196
void scaleHessian()
Scale hessian according to eqn 6.20 of Nocedal and Wright.
std::vector< double > d_gradientNew
Definition BFGSNonLinearSolver.h:167
std::vector< double > d_valueNew
Definition BFGSNonLinearSolver.h:167
std::vector< double > d_deltaXNew
Definition BFGSNonLinearSolver.h:172
bool updateSolution(const std::vector< double > &step, nonlinearSolverProblem &problem)
Update solution x -> x + step.
bool d_stepAccepted
boolean parameter for step accepteance and Wolfe conditions.
Definition BFGSNonLinearSolver.h:196
double d_trustRadiusInitial
storage for trust region parameters.
Definition BFGSNonLinearSolver.h:192
nonLinearSolver::ReturnValueType solve(nonlinearSolverProblem &problem, const std::string checkpointFileName="", const bool restart=false)
Solve non-linear problem using a modified BFGS method.
void computeNewtonStep()
Compute the Quasi-Newton Step.
std::vector< double > d_value
Definition BFGSNonLinearSolver.h:163
void load(const std::string &checkpointFileName)
Load BFGS solver state from checkpoint file.
const bool d_usePreconditioner
Definition BFGSNonLinearSolver.h:213
double d_trustRadius
Definition BFGSNonLinearSolver.h:193
double d_trustRadiusMax
Definition BFGSNonLinearSolver.h:192
void save(const std::string &checkpointFileName)
Create checkpoint file for current state of the BFGS solver.
unsigned int d_iter
storage for current bfgs iteration count
Definition BFGSNonLinearSolver.h:179
unsigned int d_numberUnknowns
Definition BFGSNonLinearSolver.h:176
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.
void computeStep()
Compute the final update step using the trust radius and whether or not the previous step was accepte...
void initializeHessian(nonlinearSolverProblem &problem)
initialize hessian, either preconditioner or identity matrix.
double d_trustRadiusMin
Definition BFGSNonLinearSolver.h:192
void updateHessian()
Update Hessian according to damped BFGS rule: Procedure 18.2 of Nocedal and Wright.
std::vector< double > d_updateVector
Definition BFGSNonLinearSolver.h:172
int d_isReset
Definition BFGSNonLinearSolver.h:210
bool d_wolfeSufficientDec
Definition BFGSNonLinearSolver.h:196
bool d_useSingleAtomSolutionsInitialGuess
Definition BFGSNonLinearSolver.h:212
bool d_wolfeSatisfied
Definition BFGSNonLinearSolver.h:197
std::vector< double > d_Srfo
Definition BFGSNonLinearSolver.h:183
bool d_hessianScaled
Definition BFGSNonLinearSolver.h:202
void computeTrustRadius(nonlinearSolverProblem &problem)
Estimate the trust radius for the next step based on the previous step and check for trust radius max...
double d_normDeltaXnew
storage for inf norm of the update step.
Definition BFGSNonLinearSolver.h:189
std::vector< double > d_hessian
storage for the hessian in current bfgs step.
Definition BFGSNonLinearSolver.h:186
MPI_Comm mpi_communicator
Definition BFGSNonLinearSolver.h:216
std::vector< double > d_gradient
Definition BFGSNonLinearSolver.h:163
void computeRFOStep()
Compute step using the Rational Function Method.
std::vector< double > d_deltaX
Definition BFGSNonLinearSolver.h:172
bool d_isCurvatureOnlyLineSearchStoppingCondition
Definition BFGSNonLinearSolver.h:205
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