Concrete class implementing Polak-Ribiere-Polyak Conjugate Gradient non-linear algebraic solver.
More...
#include <cgPRPNonLinearSolver.h>
|
| cgPRPNonLinearSolver (const unsigned int maxNumberIterations, const unsigned int debugLevel, const MPI_Comm &mpi_comm_parent, const double lineSearchTolerance=1.0e-6, const unsigned int lineSearchMaxIterations=10, const double lineSeachDampingParameter=1.0, const double maxIncrementSolLinf=1e+6, const bool isCurvatureOnlyLineSearchStoppingCondition=false) |
| Constructor.
|
|
| ~cgPRPNonLinearSolver () |
| Destructor.
|
|
nonLinearSolver::ReturnValueType | solve (nonlinearSolverProblem &problem, const std::string checkpointFileName="", const bool restart=false) |
| Solve non-linear problem using Polak-Ribiere-Polyak nonlinar conjugate gradient method.
|
|
void | save (const std::string &checkpointFileName) |
| Create checkpoint file for current state of the cg solver.
|
|
virtual | ~nonLinearSolver ()=0 |
| Destructor.
|
|
Concrete class implementing Polak-Ribiere-Polyak Conjugate Gradient non-linear algebraic solver.
- Author
- Sambit Das
◆ cgPRPNonLinearSolver()
dftfe::cgPRPNonLinearSolver::cgPRPNonLinearSolver |
( |
const unsigned int | maxNumberIterations, |
|
|
const unsigned int | debugLevel, |
|
|
const MPI_Comm & | mpi_comm_parent, |
|
|
const double | lineSearchTolerance = 1.0e-6, |
|
|
const unsigned int | lineSearchMaxIterations = 10, |
|
|
const double | lineSeachDampingParameter = 1.0, |
|
|
const double | maxIncrementSolLinf = 1e+6, |
|
|
const bool | isCurvatureOnlyLineSearchStoppingCondition = false ) |
Constructor.
- Parameters
-
maxNumberIterations | Maximum number of iterations. |
debugLevel | Debug output level: 0 - no debug output 1 - limited debug output 2 - all debug output. |
lineSearchTolerance | Tolereance required for line search convergence. |
lineSearchMaxIterations | Maximum number of iterations for the line search. |
lineSearchDampingParameter | scales the initial line search step |
◆ ~cgPRPNonLinearSolver()
dftfe::cgPRPNonLinearSolver::~cgPRPNonLinearSolver |
( |
| ) |
|
◆ computeDeltaD()
std::pair< double, double > dftfe::cgPRPNonLinearSolver::computeDeltaD |
( |
| ) |
|
|
private |
Compute delta_d and eta_p.
- Returns
- Pair containing delta_d and eta_p.
◆ computeDeltas()
void dftfe::cgPRPNonLinearSolver::computeDeltas |
( |
| ) |
|
|
private |
Compute delta new and delta mid.
◆ computeEta()
double dftfe::cgPRPNonLinearSolver::computeEta |
( |
| ) |
|
|
private |
Compute eta.
- Returns
- Value of eta.
◆ computeResidualL2Norm()
double dftfe::cgPRPNonLinearSolver::computeResidualL2Norm |
( |
| ) |
const |
|
private |
Compute residual L2-norm.
- Returns
- Value of the residual L2-norm.
◆ computeTotalNumberUnknowns()
unsigned int dftfe::cgPRPNonLinearSolver::computeTotalNumberUnknowns |
( |
| ) |
const |
|
private |
Compute the total number of unknowns in all processors.
- Returns
- Number of unknowns in all processors.
◆ initializeDirection()
void dftfe::cgPRPNonLinearSolver::initializeDirection |
( |
| ) |
|
|
private |
◆ lineSearch()
nonLinearSolver::ReturnValueType dftfe::cgPRPNonLinearSolver::lineSearch |
( |
nonlinearSolverProblem & | problem, |
|
|
const double | tolerance, |
|
|
const unsigned int | maxNumberIterations, |
|
|
const unsigned int | debugLevel, |
|
|
const std::string | checkpointFileName = "", |
|
|
const int | startingIter = -1, |
|
|
const bool | isCheckpointRestart = false ) |
|
private |
Perform line search.
- Parameters
-
problem | nonlinearSolverProblem object (functor) to compute energy and forces. |
tolerance | Tolerance (relative) required for convergence. |
maxNumberIterations | Maximum number of iterations. |
debugLevel | Debug output level: 0 - no debug output 1 - limited debug output 2 - all debug output. |
- Returns
- Return value indicating success or failure.
◆ load()
void dftfe::cgPRPNonLinearSolver::load |
( |
const std::string & | checkpointFileName | ) |
|
|
private |
Load cg solver state from checkpoint file.
◆ save()
void dftfe::cgPRPNonLinearSolver::save |
( |
const std::string & | checkpointFileName | ) |
|
|
virtual |
◆ solve()
Solve non-linear problem using Polak-Ribiere-Polyak nonlinar conjugate gradient 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.
◆ updateDirection()
void dftfe::cgPRPNonLinearSolver::updateDirection |
( |
| ) |
|
|
private |
◆ updateSolution()
bool dftfe::cgPRPNonLinearSolver::updateSolution |
( |
const double | alpha, |
|
|
const std::vector< double > & | direction, |
|
|
nonlinearSolverProblem & | problem ) |
|
private |
Update solution x -> x + \alpha direction.
- Parameters
-
alpha | Scalar value. |
direction | Direction vector. |
problem | nonlinearSolverProblem object. |
- Returns
- bool true if valid update and false if increment bound exceeded
◆ d_alphaChk
double dftfe::cgPRPNonLinearSolver::d_alphaChk |
|
private |
◆ d_beta
double dftfe::cgPRPNonLinearSolver::d_beta |
|
private |
storage for beta- the parameter for updating the conjugate direction d_beta = (d_deltaNew - d_deltaMid)/d_deltaOld
◆ d_conjugateDirection
std::vector<double> dftfe::cgPRPNonLinearSolver::d_conjugateDirection |
|
private |
storage for conjugate direction
◆ d_deltaMid
double dftfe::cgPRPNonLinearSolver::d_deltaMid |
|
private |
intermediate variable for beta computation
◆ d_deltaNew
double dftfe::cgPRPNonLinearSolver::d_deltaNew |
|
private |
intermediate variable for beta computation
◆ d_deltaOld
double dftfe::cgPRPNonLinearSolver::d_deltaOld |
|
private |
intermediate variable for beta computation
◆ d_eta
double dftfe::cgPRPNonLinearSolver::d_eta |
|
private |
◆ d_etaAlphaZeroChk
double dftfe::cgPRPNonLinearSolver::d_etaAlphaZeroChk |
|
private |
◆ d_etaChk
double dftfe::cgPRPNonLinearSolver::d_etaChk |
|
private |
◆ d_etaPChk
double dftfe::cgPRPNonLinearSolver::d_etaPChk |
|
private |
◆ d_functionalValueAfterAlphUpdateChk
double dftfe::cgPRPNonLinearSolver::d_functionalValueAfterAlphUpdateChk |
|
private |
◆ d_functionValueChk
double dftfe::cgPRPNonLinearSolver::d_functionValueChk |
|
private |
◆ d_gradient
std::vector<double> dftfe::cgPRPNonLinearSolver::d_gradient |
|
private |
storage for the gradient of the nonlinear problem in the current cg step
◆ d_gradMax
double dftfe::cgPRPNonLinearSolver::d_gradMax |
|
private |
◆ d_isCGRestartDueToLargeIncrement
bool dftfe::cgPRPNonLinearSolver::d_isCGRestartDueToLargeIncrement |
|
private |
flag which restarts the CG if large increment to the solution vector happens during line search
◆ d_isCurvatureOnlyLineSearchStoppingCondition
bool dftfe::cgPRPNonLinearSolver::d_isCurvatureOnlyLineSearchStoppingCondition |
|
private |
◆ d_iter
unsigned int dftfe::cgPRPNonLinearSolver::d_iter |
|
private |
storage for current nonlinear cg iteration count
◆ d_lineSearchDampingParameter
double dftfe::cgPRPNonLinearSolver::d_lineSearchDampingParameter |
|
private |
damping parameter (0,1] to be multiplied with the steepest descent direction, which controls the initial guess to the line search iteration.
◆ d_lineSearchMaxIterations
const unsigned int dftfe::cgPRPNonLinearSolver::d_lineSearchMaxIterations |
|
private |
maximum number of line search iterations
◆ d_lineSearchRestartIterChk
int dftfe::cgPRPNonLinearSolver::d_lineSearchRestartIterChk |
|
private |
◆ d_lineSearchTolerance
const double dftfe::cgPRPNonLinearSolver::d_lineSearchTolerance |
|
private |
line search stopping tolerance
◆ d_maxSolutionIncrementLinf
double dftfe::cgPRPNonLinearSolver::d_maxSolutionIncrementLinf |
|
private |
maximum allowed increment (measured as L_{inf}(delta x)) in solution vector beyond which CG is restarted
◆ d_numberUnknowns
unsigned int dftfe::cgPRPNonLinearSolver::d_numberUnknowns |
|
private |
storage for number of unknowns to be solved for in the nonlinear problem
◆ d_steepestDirectionOld
std::vector<double> dftfe::cgPRPNonLinearSolver::d_steepestDirectionOld |
|
private |
storage for the steepest descent direction of the nonlinear problem in the previous cg step
◆ d_unknownCountFlag
std::vector<unsigned int> dftfe::cgPRPNonLinearSolver::d_unknownCountFlag |
|
private |
Storage for vector of flags (0 or 1) with size equal to the size of the solution vector of the nonlinear problem. If the flag value is 1 for an index in the vector, the corresponding entry in the solution vector is allowed to be updated and vice-versa if flag value is 0 for an index.
◆ d_useSingleAtomSolutionsInitialGuess
bool dftfe::cgPRPNonLinearSolver::d_useSingleAtomSolutionsInitialGuess |
|
private |
◆ mpi_communicator
MPI_Comm dftfe::cgPRPNonLinearSolver::mpi_communicator |
|
private |
◆ n_mpi_processes
const unsigned int dftfe::cgPRPNonLinearSolver::n_mpi_processes |
|
private |
◆ pcout
dealii::ConditionalOStream dftfe::cgPRPNonLinearSolver::pcout |
|
private |
◆ this_mpi_process
const unsigned int dftfe::cgPRPNonLinearSolver::this_mpi_process |
|
private |
The documentation for this class was generated from the following file: