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

problem class for atomic force relaxation solver. More...

#include <geoOptIon.h>

Inheritance diagram for dftfe::geoOptIon:
dftfe::nonlinearSolverProblem

Public Member Functions

 geoOptIon (dftBase *dftPtr, const MPI_Comm &mpi_comm_parent, const bool restart=false)
 Constructor.
 
void init (const std::string &restartPath)
 initializes the data member d_relaxationFlags.
 
int run ()
 calls the atomic force relaxation solver.
 
unsigned int getNumberUnknowns () const
 Obtain number of unknowns (total number of force components to be relaxed).
 
void gradient (std::vector< double > &gradient)
 Compute function gradient (aka forces).
 
void update (const std::vector< double > &solution, const bool computeForces=true, const bool useSingleAtomSolutionsInitialGuess=false)
 Update atomic positions.
 
void save ()
 create checkpoint file for current domain bounding vectors and atomic coordinates.
 
bool isConverged () const
 check for convergence.
 
const MPI_Comm & getMPICommunicator ()
 get MPI communicator.
 
void value (std::vector< double > &functionValue)
 not implemented
 
void precondition (std::vector< double > &s, const std::vector< double > &gradient)
 not implemented
 
void solution (std::vector< double > &solution)
 not implemented
 
std::vector< unsigned int > getUnknownCountFlag () const
 not implemented
 
- Public Member Functions inherited from dftfe::nonlinearSolverProblem
 nonlinearSolverProblem ()
 Constructor.
 
virtual ~nonlinearSolverProblem ()=0
 Destructor.
 

Private Attributes

std::vector< unsigned int > d_relaxationFlags
 
std::vector< double > d_externalForceOnAtom
 
std::vector< std::vector< double > > d_atomLocationsInitial
 
std::string d_restartPath
 
std::string d_solverRestartPath
 
bool d_isRestart
 
bool d_solverRestart
 
bool d_isScfRestart
 
int d_solver
 
double d_maximumAtomForceToBeRelaxed
 maximum force component to be relaxed
 
int d_totalUpdateCalls
 total number of calls to update()
 
dftBased_dftPtr
 pointer to dft class
 
std::unique_ptr< nonLinearSolverd_nonLinearSolverPtr
 
const MPI_Comm mpi_communicator
 parallel communication objects
 
const unsigned int n_mpi_processes
 
const unsigned int this_mpi_process
 
dealii::ConditionalOStream pcout
 conditional stream object
 

Detailed Description

problem class for atomic force relaxation solver.

Author
Sambit Das

Constructor & Destructor Documentation

◆ geoOptIon()

dftfe::geoOptIon::geoOptIon ( dftBase * dftPtr,
const MPI_Comm & mpi_comm_parent,
const bool restart = false )

Constructor.

Parameters
_dftPtrpointer to dftClass
mpi_comm_parentparent mpi_communicator

Member Function Documentation

◆ getMPICommunicator()

const MPI_Comm & dftfe::geoOptIon::getMPICommunicator ( )
virtual

get MPI communicator.

Implements dftfe::nonlinearSolverProblem.

◆ getNumberUnknowns()

unsigned int dftfe::geoOptIon::getNumberUnknowns ( ) const
virtual

Obtain number of unknowns (total number of force components to be relaxed).

Returns
int Number of unknowns.

Implements dftfe::nonlinearSolverProblem.

◆ getUnknownCountFlag()

std::vector< unsigned int > dftfe::geoOptIon::getUnknownCountFlag ( ) const
virtual

not implemented

Implements dftfe::nonlinearSolverProblem.

◆ gradient()

void dftfe::geoOptIon::gradient ( std::vector< double > & gradient)
virtual

Compute function gradient (aka forces).

Parameters
gradientSTL vector for gradient values.

Implements dftfe::nonlinearSolverProblem.

◆ init()

void dftfe::geoOptIon::init ( const std::string & restartPath)

initializes the data member d_relaxationFlags.

◆ isConverged()

bool dftfe::geoOptIon::isConverged ( ) const
virtual

check for convergence.

Implements dftfe::nonlinearSolverProblem.

◆ precondition()

void dftfe::geoOptIon::precondition ( std::vector< double > & s,
const std::vector< double > & gradient )
virtual

not implemented

Implements dftfe::nonlinearSolverProblem.

◆ run()

int dftfe::geoOptIon::run ( )

calls the atomic force relaxation solver.

Currently we have option of one solver: Polak–Ribière nonlinear CG solver with secant based line search. In future releases, we will have more options like BFGS solver.

Returns
int total geometry update calls

◆ save()

void dftfe::geoOptIon::save ( )
virtual

create checkpoint file for current domain bounding vectors and atomic coordinates.

Implements dftfe::nonlinearSolverProblem.

◆ solution()

void dftfe::geoOptIon::solution ( std::vector< double > & solution)
virtual

not implemented

Implements dftfe::nonlinearSolverProblem.

◆ update()

void dftfe::geoOptIon::update ( const std::vector< double > & solution,
const bool computeForces = true,
const bool useSingleAtomSolutionsInitialGuess = false )
virtual

Update atomic positions.

Parameters
solutiondisplacement of the atoms with respect to their current position. The size of the solution vector is equal to the number of unknowns.

Implements dftfe::nonlinearSolverProblem.

◆ value()

void dftfe::geoOptIon::value ( std::vector< double > & functionValue)
virtual

not implemented

Implements dftfe::nonlinearSolverProblem.

Member Data Documentation

◆ d_atomLocationsInitial

std::vector<std::vector<double> > dftfe::geoOptIon::d_atomLocationsInitial
private

◆ d_dftPtr

dftBase* dftfe::geoOptIon::d_dftPtr
private

pointer to dft class

◆ d_externalForceOnAtom

std::vector<double> dftfe::geoOptIon::d_externalForceOnAtom
private

◆ d_isRestart

bool dftfe::geoOptIon::d_isRestart
private

◆ d_isScfRestart

bool dftfe::geoOptIon::d_isScfRestart
private

◆ d_maximumAtomForceToBeRelaxed

double dftfe::geoOptIon::d_maximumAtomForceToBeRelaxed
private

maximum force component to be relaxed

◆ d_nonLinearSolverPtr

std::unique_ptr<nonLinearSolver> dftfe::geoOptIon::d_nonLinearSolverPtr
private

◆ d_relaxationFlags

std::vector<unsigned int> dftfe::geoOptIon::d_relaxationFlags
private

storage for relaxation flags and external force components for each global atom. each atom has three flags corresponding to three components (0- no relax, 1- relax) and three external force components

◆ d_restartPath

std::string dftfe::geoOptIon::d_restartPath
private

◆ d_solver

int dftfe::geoOptIon::d_solver
private

◆ d_solverRestart

bool dftfe::geoOptIon::d_solverRestart
private

◆ d_solverRestartPath

std::string dftfe::geoOptIon::d_solverRestartPath
private

◆ d_totalUpdateCalls

int dftfe::geoOptIon::d_totalUpdateCalls
private

total number of calls to update()

◆ mpi_communicator

const MPI_Comm dftfe::geoOptIon::mpi_communicator
private

parallel communication objects

◆ n_mpi_processes

const unsigned int dftfe::geoOptIon::n_mpi_processes
private

◆ pcout

dealii::ConditionalOStream dftfe::geoOptIon::pcout
private

conditional stream object

◆ this_mpi_process

const unsigned int dftfe::geoOptIon::this_mpi_process
private

The documentation for this class was generated from the following file: