|
DFT-FE 1.3.0-pre
Density Functional Theory With Finite-Elements
|
problem class for atomic force relaxation solver. More...
#include <geoOptIon.h>
Public Member Functions | |
| geoOptIon (dftBase *dftPtr, const MPI_Comm &mpi_comm_parent, const bool restart=false) | |
| Constructor. | |
| void | init (const std::string &restartPath, const dftfe::Int cycleId) |
| initializes the data member d_relaxationFlags. | |
| dftfe::Int | run () |
| calls the atomic force relaxation solver. | |
| dftfe::uInt | 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< dftfe::uInt > | getUnknownCountFlag () const |
| not implemented | |
Public Member Functions inherited from dftfe::nonlinearSolverProblem | |
| nonlinearSolverProblem () | |
| Constructor. | |
| virtual | ~nonlinearSolverProblem ()=0 |
| Destructor. | |
Private Attributes | |
| std::vector< dftfe::uInt > | 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 |
| dftfe::Int | d_solver |
| double | d_maximumAtomForceToBeRelaxed |
| maximum force component to be relaxed | |
| dftfe::Int | d_totalUpdateCalls |
| total number of calls to update() | |
| dftfe::Int | d_cycle |
| dftBase * | d_dftPtr |
| pointer to dft class | |
| std::unique_ptr< nonLinearSolver > | d_nonLinearSolverPtr |
| const MPI_Comm | mpi_communicator |
| parallel communication objects | |
| const dftfe::uInt | n_mpi_processes |
| const dftfe::uInt | this_mpi_process |
| dealii::ConditionalOStream | pcout |
| conditional stream object | |
problem class for atomic force relaxation solver.
| dftfe::geoOptIon::geoOptIon | ( | dftBase * | dftPtr, |
| const MPI_Comm & | mpi_comm_parent, | ||
| const bool | restart = false ) |
Constructor.
| _dftPtr | pointer to dftClass |
| mpi_comm_parent | parent mpi_communicator |
|
virtual |
get MPI communicator.
Implements dftfe::nonlinearSolverProblem.
|
virtual |
Obtain number of unknowns (total number of force components to be relaxed).
Implements dftfe::nonlinearSolverProblem.
|
virtual |
not implemented
Implements dftfe::nonlinearSolverProblem.
|
virtual |
Compute function gradient (aka forces).
| gradient | STL vector for gradient values. |
Implements dftfe::nonlinearSolverProblem.
| void dftfe::geoOptIon::init | ( | const std::string & | restartPath, |
| const dftfe::Int | cycleId ) |
initializes the data member d_relaxationFlags.
|
virtual |
check for convergence.
Implements dftfe::nonlinearSolverProblem.
|
virtual |
not implemented
Implements dftfe::nonlinearSolverProblem.
| dftfe::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.
|
virtual |
create checkpoint file for current domain bounding vectors and atomic coordinates.
Implements dftfe::nonlinearSolverProblem.
|
virtual |
not implemented
Implements dftfe::nonlinearSolverProblem.
|
virtual |
Update atomic positions.
| solution | displacement 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.
|
virtual |
not implemented
Implements dftfe::nonlinearSolverProblem.
|
private |
|
private |
|
private |
pointer to dft class
|
private |
|
private |
|
private |
|
private |
maximum force component to be relaxed
|
private |
|
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
|
private |
|
private |
|
private |
|
private |
|
private |
total number of calls to update()
|
private |
parallel communication objects
|
private |
|
private |
conditional stream object
|
private |