DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftfe::cgPRPNonLinearSolver Class Reference

Concrete class implementing Polak-Ribiere-Polyak Conjugate Gradient non-linear algebraic solver. More...

#include <cgPRPNonLinearSolver.h>

Inheritance diagram for dftfe::cgPRPNonLinearSolver:
dftfe::nonLinearSolver

Public Member Functions

 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.
 
- Public Member Functions inherited from dftfe::nonLinearSolver
virtual ~nonLinearSolver ()=0
 Destructor.
 

Private Member Functions

void initializeDirection ()
 Initialize direction.
 
nonLinearSolver::ReturnValueType 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)
 Perform line search.
 
std::pair< double, double > computeDeltaD ()
 Compute delta_d and eta_p.
 
double computeEta ()
 Compute eta.
 
void computeDeltas ()
 Compute delta new and delta mid.
 
void updateDirection ()
 Update direction.
 
double computeResidualL2Norm () const
 Compute residual L2-norm.
 
unsigned int computeTotalNumberUnknowns () const
 Compute the total number of unknowns in all processors.
 
bool updateSolution (const double alpha, const std::vector< double > &direction, nonlinearSolverProblem &problem)
 Update solution x -> x + \alpha direction.
 
void load (const std::string &checkpointFileName)
 Load cg solver state from checkpoint file.
 

Private Attributes

std::vector< double > d_conjugateDirection
 storage for conjugate direction
 
std::vector< double > d_gradient
 storage for the gradient of the nonlinear problem in the current cg step
 
std::vector< double > d_steepestDirectionOld
 
double d_deltaNew
 intermediate variable for beta computation
 
double d_deltaMid
 intermediate variable for beta computation
 
double d_deltaOld
 intermediate variable for beta computation
 
double d_beta
 
double d_gradMax
 
unsigned int d_numberUnknowns
 storage for number of unknowns to be solved for in the nonlinear problem
 
unsigned int d_iter
 storage for current nonlinear cg iteration count
 
std::vector< unsigned int > d_unknownCountFlag
 
const double d_lineSearchTolerance
 line search stopping tolerance
 
const unsigned int d_lineSearchMaxIterations
 maximum number of line search iterations
 
double d_lineSearchDampingParameter
 
bool d_isCGRestartDueToLargeIncrement
 
double d_maxSolutionIncrementLinf
 
double d_alphaChk
 line search data
 
double d_etaPChk
 line search data
 
double d_etaChk
 line search data
 
double d_eta
 line search data
 
double d_etaAlphaZeroChk
 line search data
 
double d_functionValueChk
 line search data
 
double d_functionalValueAfterAlphUpdateChk
 line search data
 
int d_lineSearchRestartIterChk
 line search iter
 
bool d_useSingleAtomSolutionsInitialGuess
 
bool d_isCurvatureOnlyLineSearchStoppingCondition
 
MPI_Comm mpi_communicator
 
const unsigned int n_mpi_processes
 
const unsigned int this_mpi_process
 
dealii::ConditionalOStream pcout
 

Additional Inherited Members

- Public Types inherited from dftfe::nonLinearSolver
enum  ReturnValueType {
  SUCCESS = 0 , FAILURE , LINESEARCH_FAILED , MAX_ITER_REACHED ,
  RESTART
}
 
- Protected Member Functions inherited from dftfe::nonLinearSolver
 nonLinearSolver (const unsigned int debugLevel, const unsigned int maxNumberIterations, const double tolerance=0)
 Constructor.
 
double getTolerance () const
 Get tolerance.
 
unsigned int getMaximumNumberIterations () const
 Get maximum number of iterations.
 
unsigned int getDebugLevel () const
 Get debug level.
 
- Protected Attributes inherited from dftfe::nonLinearSolver
const unsigned int d_debugLevel
 controls the verbosity of the printing
 
const unsigned int d_maxNumberIterations
 maximum number of nonlinear solve iterations
 
const double d_tolerance
 nonlinear solve stopping tolerance
 

Detailed Description

Concrete class implementing Polak-Ribiere-Polyak Conjugate Gradient non-linear algebraic solver.

Author
Sambit Das

Constructor & Destructor Documentation

◆ 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
maxNumberIterationsMaximum number of iterations.
debugLevelDebug output level: 0 - no debug output 1 - limited debug output 2 - all debug output.
lineSearchToleranceTolereance required for line search convergence.
lineSearchMaxIterationsMaximum number of iterations for the line search.
lineSearchDampingParameterscales the initial line search step

◆ ~cgPRPNonLinearSolver()

dftfe::cgPRPNonLinearSolver::~cgPRPNonLinearSolver ( )

Destructor.

Member Function Documentation

◆ 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

Initialize direction.

◆ 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
problemnonlinearSolverProblem object (functor) to compute energy and forces.
toleranceTolerance (relative) required for convergence.
maxNumberIterationsMaximum number of iterations.
debugLevelDebug 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

Create checkpoint file for current state of the cg solver.

Implements dftfe::nonLinearSolver.

◆ solve()

nonLinearSolver::ReturnValueType dftfe::cgPRPNonLinearSolver::solve ( nonlinearSolverProblem & problem,
const std::string checkpointFileName = "",
const bool restart = false )
virtual

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

Update direction.

◆ updateSolution()

bool dftfe::cgPRPNonLinearSolver::updateSolution ( const double alpha,
const std::vector< double > & direction,
nonlinearSolverProblem & problem )
private

Update solution x -> x + \alpha direction.

Parameters
alphaScalar value.
directionDirection vector.
problemnonlinearSolverProblem object.
Returns
bool true if valid update and false if increment bound exceeded

Member Data Documentation

◆ d_alphaChk

double dftfe::cgPRPNonLinearSolver::d_alphaChk
private

line search data

◆ 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

line search data

◆ d_etaAlphaZeroChk

double dftfe::cgPRPNonLinearSolver::d_etaAlphaZeroChk
private

line search data

◆ d_etaChk

double dftfe::cgPRPNonLinearSolver::d_etaChk
private

line search data

◆ d_etaPChk

double dftfe::cgPRPNonLinearSolver::d_etaPChk
private

line search data

◆ d_functionalValueAfterAlphUpdateChk

double dftfe::cgPRPNonLinearSolver::d_functionalValueAfterAlphUpdateChk
private

line search data

◆ d_functionValueChk

double dftfe::cgPRPNonLinearSolver::d_functionValueChk
private

line search data

◆ 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

line search iter

◆ 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: