DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
|
wrapper class for dftfe More...
#include <dftfeWrapper.h>
Public Member Functions | |
dftfeWrapper () | |
empty constructor | |
dftfeWrapper (const std::string parameter_file, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false) | |
constructor based on input parameter_file | |
dftfeWrapper (const std::string parameter_file, const std::string restartCoordsFile, const std::string restartDomainVectorsFile, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false, const bool isScfRestart=true) | |
constructor based on input parameter_file and restart coordinates and domain vectors file paths | |
dftfeWrapper (const MPI_Comm &mpi_comm_parent, const bool useDevice, const std::vector< std::vector< double > > atomicPositionsCart, const std::vector< unsigned int > atomicNumbers, const std::vector< std::vector< double > > cell, const std::vector< bool > pbc, const std::vector< unsigned int > mpGrid=std::vector< unsigned int >{1, 1, 1}, const std::vector< bool > mpGridShift=std::vector< bool >{false, false, false}, const bool spinPolarizedDFT=false, const double startMagnetization=0.0, const double fermiDiracSmearingTemp=500.0, const unsigned int npkpt=0, const double meshSize=0.8, const double scfMixingParameter=0.2, const int verbosity=-1, const bool setDeviceToMPITaskBindingInternally=false) | |
constructor based on input list of atomic coordinates, list of atomic numbers,cell, boundary conditions, Monkhorst-Pack k-point grid, and other optional parameters. This constructor currently only sets up GGA PBE pseudopotential DFT calculations using ONCV pseudopotentials in .upf format (read from DFTFE_PSP_PATH folder provided as an environment variable). The pseudpotential directory must contain files in the format: AtomicSymbol.upf | |
~dftfeWrapper () | |
void | reinit (const std::string parameter_file, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false) |
clear and reinitialize based on input parameter_file | |
void | reinit (const std::string parameter_file, const std::string restartCoordsFile, const std::string restartDomainVectorsFile, const MPI_Comm &mpi_comm_parent, const bool printParams=false, const bool setDeviceToMPITaskBindingInternally=false, const std::string solverMode="GS", const std::string restartFilesPath=".", const int _verbosity=1, const bool useDevice=false, const bool isScfRestart=true) |
clear and reinitialize based on input parameter_file and restart coordinates and domain vectors file paths | |
void | reinit (const MPI_Comm &mpi_comm_parent, const bool useDevice, const std::vector< std::vector< double > > atomicPositionsCart, const std::vector< unsigned int > atomicNumbers, const std::vector< std::vector< double > > cell, const std::vector< bool > pbc, const std::vector< unsigned int > mpGrid=std::vector< unsigned int >{1, 1, 1}, const std::vector< bool > mpGridShift=std::vector< bool >{false, false, false}, const bool spinPolarizedDFT=false, const double startMagnetization=0.0, const double fermiDiracSmearingTemp=500.0, const unsigned int npkpt=0, const double meshSize=0.8, const double scfMixingParameter=0.2, const int verbosity=-1, const bool setDeviceToMPITaskBindingInternally=false) |
void | clear () |
void | run () |
Legacy function (to be deprecated) | |
void | writeMesh () |
Calls dftBasepointer public function writeMesh(). Here the inital density and mesh are stored in a file. Useful for visluiing meshes without running solve. | |
std::tuple< double, bool, double > | computeDFTFreeEnergy (const bool computeIonForces=true, const bool computeCellStress=false) |
solve ground-state and return DFT free energy which is sum of internal energy and negative of electronic entropic energy (in Hartree units) | |
void | computeStress () |
double | getDFTFreeEnergy () const |
Get DFT free energy (in Hartree units). This function can only be called after calling computeDFTFreeEnergy. | |
double | getElectronicEntropicEnergy () const |
Get electronic entropic energy (in Hartree units). This function can only be called after calling computeDFTFreeEnergy. | |
std::vector< std::vector< double > > | getForcesAtoms () const |
Get ionic forces: negative of gradient of DFT free energy with respect to ionic positions (in Hartree/Bohr units). This function can only be called after calling computeDFTFreeEnergy. | |
std::vector< std::vector< double > > | getCellStress () const |
Get cell stress: negative of gradient of DFT free energy with respect to affine strain components scaled by volume (Hartree/Bohr^3) units. This function can only be called after calling computeDFTFreeEnergy. | |
void | updateAtomPositions (const std::vector< std::vector< double > > atomsDisplacements) |
update atom positions and reinitialize all related data-structures | |
void | deformCell (const std::vector< std::vector< double > > deformationGradient) |
Deforms the cell by applying the given affine deformation gradient and reinitializes the underlying data-structures. | |
std::vector< std::vector< double > > | getAtomPositionsCart () const |
Gets the current atom Positions in cartesian form (in Bohr units) (origin at corner of cell against which the cell vectors are defined) | |
std::vector< std::vector< double > > | getAtomPositionsFrac () const |
Gets the current atom Positions in fractional form (only applicable for periodic and semi-periodic BCs). CAUTION: during relaxation and MD fractional coordinates may have negaive values. | |
std::vector< std::vector< double > > | getCell () const |
Gets the current cell vectors. | |
std::vector< bool > | getPBC () const |
Gets the boundary conditions for each cell vector direction. | |
std::vector< int > | getAtomicNumbers () const |
Gets the atomic numbers vector. | |
std::vector< int > | getValenceElectronNumbers () const |
Gets the number of valence electrons for each atom. | |
void | writeDomainAndAtomCoordinates (const std::string Path) const |
writes the current domain bounding vectors and atom coordinates to files for structural optimization and dynamics restarts. The coordinates are stored as: 1. fractional for semi-periodic/periodic 2. Cartesian for non-periodic. | |
dftBase * | getDftfeBasePtr () |
dftParameters * | getDftfeParamsPtr () |
Static Public Member Functions | |
static void | globalHandlesInitialize (const MPI_Comm &mpi_comm_world) |
must be called only once at start of program from all processors after calling MPI_Init | |
static void | globalHandlesFinalize () |
must be called only once at end of program from all processors but before calling MPI_Finalize | |
Private Member Functions | |
void | createScratchFolder () |
void | initialize (const bool setDeviceToMPITaskBindingInternally, const bool useDevice) |
Private Attributes | |
MPI_Comm | d_mpi_comm_parent |
dftBase * | d_dftfeBasePtr |
dftParameters * | d_dftfeParamsPtr |
std::string | d_scratchFolderName |
bool | d_isDeviceToMPITaskBindingSetInternally |
wrapper class for dftfe
dftfe::dftfeWrapper::dftfeWrapper | ( | ) |
empty constructor
dftfe::dftfeWrapper::dftfeWrapper | ( | const std::string | parameter_file, |
const MPI_Comm & | mpi_comm_parent, | ||
const bool | printParams = false, | ||
const bool | setDeviceToMPITaskBindingInternally = false, | ||
const std::string | solverMode = "GS", | ||
const std::string | restartFilesPath = ".", | ||
const int | _verbosity = 1, | ||
const bool | useDevice = false ) |
constructor based on input parameter_file
dftfe::dftfeWrapper::dftfeWrapper | ( | const std::string | parameter_file, |
const std::string | restartCoordsFile, | ||
const std::string | restartDomainVectorsFile, | ||
const MPI_Comm & | mpi_comm_parent, | ||
const bool | printParams = false, | ||
const bool | setDeviceToMPITaskBindingInternally = false, | ||
const std::string | solverMode = "GS", | ||
const std::string | restartFilesPath = ".", | ||
const int | _verbosity = 1, | ||
const bool | useDevice = false, | ||
const bool | isScfRestart = true ) |
constructor based on input parameter_file and restart coordinates and domain vectors file paths
dftfe::dftfeWrapper::dftfeWrapper | ( | const MPI_Comm & | mpi_comm_parent, |
const bool | useDevice, | ||
const std::vector< std::vector< double > > | atomicPositionsCart, | ||
const std::vector< unsigned int > | atomicNumbers, | ||
const std::vector< std::vector< double > > | cell, | ||
const std::vector< bool > | pbc, | ||
const std::vector< unsigned int > | mpGrid = std::vector< unsigned int >{1, 1, 1}, | ||
const std::vector< bool > | mpGridShift = std::vector< bool >{false, false, false}, | ||
const bool | spinPolarizedDFT = false, | ||
const double | startMagnetization = 0.0, | ||
const double | fermiDiracSmearingTemp = 500.0, | ||
const unsigned int | npkpt = 0, | ||
const double | meshSize = 0.8, | ||
const double | scfMixingParameter = 0.2, | ||
const int | verbosity = -1, | ||
const bool | setDeviceToMPITaskBindingInternally = false ) |
constructor based on input list of atomic coordinates, list of atomic numbers,cell, boundary conditions, Monkhorst-Pack k-point grid, and other optional parameters. This constructor currently only sets up GGA PBE pseudopotential DFT calculations using ONCV pseudopotentials in .upf format (read from DFTFE_PSP_PATH folder provided as an environment variable). The pseudpotential directory must contain files in the format: AtomicSymbol.upf
[in] | mpi_comm_parent | mpi communicator to be used by the dftfeWrapper. |
[in] | useDevice | toggle use of Device accelerated DFT-FE |
[in] | atomicPositionsCart | vector of atomic positions for each atom (in Bohr units), Origin is at cell corner |
[in] | atomicNumbers | vector of atomic numbers |
[in] | cell | 3 \times 3 matrix in Bohr units, cell[i] denotes the ith cell vector. DFT-FE requires the cell vectors to form a right-handed coordinate system i.e. dotProduct(crossProduct(cell[0],cell[1]),cell[2])>0 |
[in] | pbc | vector of bools denoting periodic boundary conditions along the three cell vectors, false denotes non-periodic and true is periodic |
[in] | mpgrid | vector of Monkhorst-Pack grid points along the reciprocal lattice vector directions for sampling the Brillouin zone along periodic directions. Default value is a Gamma point. |
[in] | mpgridShift | vector of bools where false denotes no shift and true denotes shift by half the Monkhost-Pack grid spacing. Default value is no shift. |
[in] | spinPolarizedDFT | toggles spin-polarized DFT calculations. Default value is false |
[in] | startMagnetization | Starting magnetization to be used for spin-polarized DFT calculations (must be between -0.5 and +0.5). Corresponding magnetization per simulation domain will be (2 x START MAGNETIZATION x Number of electrons) in Bohr magneton units. |
[in] | fermiDiracSmearingTemp | Fermi-Dirac smearing temperature in Kelvin. Default value is 500 K. |
[in] | npkpt | Number of groups of MPI tasks across which the work load of the irreducible k-points is parallelised. npkpt must be a divisor of total number of MPI tasks. Default value of 0 internally sets npkt to an heuristically determined value. |
[in] | meshSize | Finite-element mesh size around the atoms in Bohr units. The default value of 0.8 is sufficient to achieve chemical accuracy in energy (0.1 mHa/atom discretization error) and forces (0.1 mHa/Bohr discretization error) for the ONCV pseudo-dojo pseudopotentials. Note that this function assumes a sixth order finite-element interpolating polynomial |
[in] | scfMixingParameter | mixing paramter for SCF fixed point iteration. Currently the Anderson mixing strategy is used. |
[in] | verbosity | printing verbosity. Default value is -1: no printing |
[in] | setDeviceToMPITaskBindingInternally | This option is only valid for Device runs. If set to true Device to MPI task binding is set inside the DFT-FE code. Default behaviour is false which assumes the binding has been externally set. |
dftfe::dftfeWrapper::~dftfeWrapper | ( | ) |
void dftfe::dftfeWrapper::clear | ( | ) |
std::tuple< double, bool, double > dftfe::dftfeWrapper::computeDFTFreeEnergy | ( | const bool | computeIonForces = true, |
const bool | computeCellStress = false ) |
solve ground-state and return DFT free energy which is sum of internal energy and negative of electronic entropic energy (in Hartree units)
void dftfe::dftfeWrapper::computeStress | ( | ) |
|
private |
void dftfe::dftfeWrapper::deformCell | ( | const std::vector< std::vector< double > > | deformationGradient | ) |
Deforms the cell by applying the given affine deformation gradient and reinitializes the underlying data-structures.
[in] | deformationGradient | deformation gradient matrix given by F[i][j]=\frac{\partial x_i}{\partial X_j} |
std::vector< int > dftfe::dftfeWrapper::getAtomicNumbers | ( | ) | const |
Gets the atomic numbers vector.
std::vector< std::vector< double > > dftfe::dftfeWrapper::getAtomPositionsCart | ( | ) | const |
Gets the current atom Positions in cartesian form (in Bohr units) (origin at corner of cell against which the cell vectors are defined)
std::vector< std::vector< double > > dftfe::dftfeWrapper::getAtomPositionsFrac | ( | ) | const |
Gets the current atom Positions in fractional form (only applicable for periodic and semi-periodic BCs). CAUTION: during relaxation and MD fractional coordinates may have negaive values.
std::vector< std::vector< double > > dftfe::dftfeWrapper::getCell | ( | ) | const |
Gets the current cell vectors.
std::vector< std::vector< double > > dftfe::dftfeWrapper::getCellStress | ( | ) | const |
Get cell stress: negative of gradient of DFT free energy with respect to affine strain components scaled by volume (Hartree/Bohr^3) units. This function can only be called after calling computeDFTFreeEnergy.
dftBase * dftfe::dftfeWrapper::getDftfeBasePtr | ( | ) |
dftParameters * dftfe::dftfeWrapper::getDftfeParamsPtr | ( | ) |
double dftfe::dftfeWrapper::getDFTFreeEnergy | ( | ) | const |
Get DFT free energy (in Hartree units). This function can only be called after calling computeDFTFreeEnergy.
double dftfe::dftfeWrapper::getElectronicEntropicEnergy | ( | ) | const |
Get electronic entropic energy (in Hartree units). This function can only be called after calling computeDFTFreeEnergy.
std::vector< std::vector< double > > dftfe::dftfeWrapper::getForcesAtoms | ( | ) | const |
Get ionic forces: negative of gradient of DFT free energy with respect to ionic positions (in Hartree/Bohr units). This function can only be called after calling computeDFTFreeEnergy.
std::vector< bool > dftfe::dftfeWrapper::getPBC | ( | ) | const |
Gets the boundary conditions for each cell vector direction.
std::vector< int > dftfe::dftfeWrapper::getValenceElectronNumbers | ( | ) | const |
Gets the number of valence electrons for each atom.
|
static |
must be called only once at end of program from all processors but before calling MPI_Finalize
|
static |
must be called only once at start of program from all processors after calling MPI_Init
|
private |
void dftfe::dftfeWrapper::reinit | ( | const MPI_Comm & | mpi_comm_parent, |
const bool | useDevice, | ||
const std::vector< std::vector< double > > | atomicPositionsCart, | ||
const std::vector< unsigned int > | atomicNumbers, | ||
const std::vector< std::vector< double > > | cell, | ||
const std::vector< bool > | pbc, | ||
const std::vector< unsigned int > | mpGrid = std::vector< unsigned int >{1, 1, 1}, | ||
const std::vector< bool > | mpGridShift = std::vector< bool >{false, false, false}, | ||
const bool | spinPolarizedDFT = false, | ||
const double | startMagnetization = 0.0, | ||
const double | fermiDiracSmearingTemp = 500.0, | ||
const unsigned int | npkpt = 0, | ||
const double | meshSize = 0.8, | ||
const double | scfMixingParameter = 0.2, | ||
const int | verbosity = -1, | ||
const bool | setDeviceToMPITaskBindingInternally = false ) |
void dftfe::dftfeWrapper::reinit | ( | const std::string | parameter_file, |
const MPI_Comm & | mpi_comm_parent, | ||
const bool | printParams = false, | ||
const bool | setDeviceToMPITaskBindingInternally = false, | ||
const std::string | solverMode = "GS", | ||
const std::string | restartFilesPath = ".", | ||
const int | _verbosity = 1, | ||
const bool | useDevice = false ) |
clear and reinitialize based on input parameter_file
void dftfe::dftfeWrapper::reinit | ( | const std::string | parameter_file, |
const std::string | restartCoordsFile, | ||
const std::string | restartDomainVectorsFile, | ||
const MPI_Comm & | mpi_comm_parent, | ||
const bool | printParams = false, | ||
const bool | setDeviceToMPITaskBindingInternally = false, | ||
const std::string | solverMode = "GS", | ||
const std::string | restartFilesPath = ".", | ||
const int | _verbosity = 1, | ||
const bool | useDevice = false, | ||
const bool | isScfRestart = true ) |
clear and reinitialize based on input parameter_file and restart coordinates and domain vectors file paths
void dftfe::dftfeWrapper::run | ( | ) |
Legacy function (to be deprecated)
void dftfe::dftfeWrapper::updateAtomPositions | ( | const std::vector< std::vector< double > > | atomsDisplacements | ) |
update atom positions and reinitialize all related data-structures
[in] | atomsDisplacements | vector of displacements for each atom (in Bohr units) |
void dftfe::dftfeWrapper::writeDomainAndAtomCoordinates | ( | const std::string | Path | ) | const |
writes the current domain bounding vectors and atom coordinates to files for structural optimization and dynamics restarts. The coordinates are stored as: 1. fractional for semi-periodic/periodic 2. Cartesian for non-periodic.
[in] | Path | The folder path to store the atom coordinates required during restart. |
void dftfe::dftfeWrapper::writeMesh | ( | ) |
Calls dftBasepointer public function writeMesh(). Here the inital density and mesh are stored in a file. Useful for visluiing meshes without running solve.
|
private |
|
private |
|
private |
|
private |
|
private |