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

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.
 
dftBasegetDftfeBasePtr ()
 
dftParametersgetDftfeParamsPtr ()
 

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
 
dftBased_dftfeBasePtr
 
dftParametersd_dftfeParamsPtr
 
std::string d_scratchFolderName
 
bool d_isDeviceToMPITaskBindingSetInternally
 

Detailed Description

wrapper class for dftfe

Author
Sambit Das

Constructor & Destructor Documentation

◆ dftfeWrapper() [1/4]

dftfe::dftfeWrapper::dftfeWrapper ( )

empty constructor

◆ dftfeWrapper() [2/4]

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

◆ dftfeWrapper() [3/4]

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

◆ dftfeWrapper() [4/4]

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

Parameters
[in]mpi_comm_parentmpi communicator to be used by the dftfeWrapper.
[in]useDevicetoggle use of Device accelerated DFT-FE
[in]atomicPositionsCartvector of atomic positions for each atom (in Bohr units), Origin is at cell corner
[in]atomicNumbersvector of atomic numbers
[in]cell3 \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]pbcvector of bools denoting periodic boundary conditions along the three cell vectors, false denotes non-periodic and true is periodic
[in]mpgridvector 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]mpgridShiftvector of bools where false denotes no shift and true denotes shift by half the Monkhost-Pack grid spacing. Default value is no shift.
[in]spinPolarizedDFTtoggles spin-polarized DFT calculations. Default value is false
[in]startMagnetizationStarting 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]fermiDiracSmearingTempFermi-Dirac smearing temperature in Kelvin. Default value is 500 K.
[in]npkptNumber 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]meshSizeFinite-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]scfMixingParametermixing paramter for SCF fixed point iteration. Currently the Anderson mixing strategy is used.
[in]verbosityprinting verbosity. Default value is -1: no printing
[in]setDeviceToMPITaskBindingInternallyThis 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.

◆ ~dftfeWrapper()

dftfe::dftfeWrapper::~dftfeWrapper ( )

Member Function Documentation

◆ clear()

void dftfe::dftfeWrapper::clear ( )

◆ computeDFTFreeEnergy()

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)

Returns
tuple of ground-state energy, boolean flag on whether scf converged, and L2 norm of residual electron-density of the last SCF iteration

◆ computeStress()

void dftfe::dftfeWrapper::computeStress ( )

◆ createScratchFolder()

void dftfe::dftfeWrapper::createScratchFolder ( )
private

◆ deformCell()

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.

Parameters
[in]deformationGradientdeformation gradient matrix given by F[i][j]=\frac{\partial x_i}{\partial X_j}

◆ getAtomicNumbers()

std::vector< int > dftfe::dftfeWrapper::getAtomicNumbers ( ) const

Gets the atomic numbers vector.

Returns
vector of atomic numbers

◆ getAtomPositionsCart()

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)

Returns
array of coords for each atom

◆ getAtomPositionsFrac()

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.

Returns
array of coords for each atom

◆ getCell()

std::vector< std::vector< double > > dftfe::dftfeWrapper::getCell ( ) const

Gets the current cell vectors.

Returns
3 \times 3 matrix, cell[i][j] corresponds to jth component of ith cell vector (in Bohr units)

◆ getCellStress()

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.

Returns
cell stress 3 \times 3 matrix given by sigma[i][j]=\frac{1}{\Omega}\frac{\partial E}{\partial \epsilon_{ij}}

◆ getDftfeBasePtr()

dftBase * dftfe::dftfeWrapper::getDftfeBasePtr ( )

◆ getDftfeParamsPtr()

dftParameters * dftfe::dftfeWrapper::getDftfeParamsPtr ( )

◆ getDFTFreeEnergy()

double dftfe::dftfeWrapper::getDFTFreeEnergy ( ) const

Get DFT free energy (in Hartree units). This function can only be called after calling computeDFTFreeEnergy.

◆ getElectronicEntropicEnergy()

double dftfe::dftfeWrapper::getElectronicEntropicEnergy ( ) const

Get electronic entropic energy (in Hartree units). This function can only be called after calling computeDFTFreeEnergy.

◆ getForcesAtoms()

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.

Returns
vector of forces on each atom

◆ getPBC()

std::vector< bool > dftfe::dftfeWrapper::getPBC ( ) const

Gets the boundary conditions for each cell vector direction.

Returns
vector of bools, false denotes non-periodic BC and true denotes periodic BC

◆ getValenceElectronNumbers()

std::vector< int > dftfe::dftfeWrapper::getValenceElectronNumbers ( ) const

Gets the number of valence electrons for each atom.

Returns
array of number of valence for each atom

◆ globalHandlesFinalize()

static void dftfe::dftfeWrapper::globalHandlesFinalize ( )
static

must be called only once at end of program from all processors but before calling MPI_Finalize

◆ globalHandlesInitialize()

static void dftfe::dftfeWrapper::globalHandlesInitialize ( const MPI_Comm & mpi_comm_world)
static

must be called only once at start of program from all processors after calling MPI_Init

◆ initialize()

void dftfe::dftfeWrapper::initialize ( const bool setDeviceToMPITaskBindingInternally,
const bool useDevice )
private

◆ reinit() [1/3]

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 )

◆ reinit() [2/3]

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

◆ reinit() [3/3]

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

◆ run()

void dftfe::dftfeWrapper::run ( )

Legacy function (to be deprecated)

◆ updateAtomPositions()

void dftfe::dftfeWrapper::updateAtomPositions ( const std::vector< std::vector< double > > atomsDisplacements)

update atom positions and reinitialize all related data-structures

Parameters
[in]atomsDisplacementsvector of displacements for each atom (in Bohr units)

◆ writeDomainAndAtomCoordinates()

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.

Parameters
[in]PathThe folder path to store the atom coordinates required during restart.

◆ writeMesh()

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.

Member Data Documentation

◆ d_dftfeBasePtr

dftBase* dftfe::dftfeWrapper::d_dftfeBasePtr
private

◆ d_dftfeParamsPtr

dftParameters* dftfe::dftfeWrapper::d_dftfeParamsPtr
private

◆ d_isDeviceToMPITaskBindingSetInternally

bool dftfe::dftfeWrapper::d_isDeviceToMPITaskBindingSetInternally
private

◆ d_mpi_comm_parent

MPI_Comm dftfe::dftfeWrapper::d_mpi_comm_parent
private

◆ d_scratchFolderName

std::string dftfe::dftfeWrapper::d_scratchFolderName
private

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