DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
|
#include <AtomCenteredSphericalFunctionContainer.h>
Public Member Functions | |
void | init (const std::vector< unsigned int > &atomicNumbers, const std::map< std::pair< unsigned int, unsigned int >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > &listOfSphericalFunctions) |
Initialises the class with the atomicNumbers of various atoms and the AtomCenteredSphericalFn of various spherical functions. This function is only called once per run. | |
void | initaliseCoordinates (const std::vector< double > &atomCoords, const std::vector< std::vector< double > > &periodicCoords, const std::vector< int > &imageIds) |
Initialises the position of atoms, the image posisiton and image ids after every update of atom positions. | |
unsigned int | getNumAtomCentersSize () |
Returns the number of atoms present in domain. | |
const std::vector< double > & | getAtomCoordinates () const |
Returns the cooridnates of atom present in domain. | |
const std::map< unsigned int, std::vector< double > > & | getPeriodicImageCoordinatesList () const |
Returns the map of atomId vs vector of image coordinates. | |
unsigned int | getTotalNumberOfSphericalFunctionsPerAtom (unsigned int atomicNumber) |
Returns the he number of total spherical functions indexed by {ilm} associated with an atomic number. If the atomic number does not exist, it returns a zero. | |
unsigned int | getTotalNumberOfRadialSphericalFunctionsPerAtom (unsigned int atomicNumber) |
Returns the he number of radial spherical functions indexed by {i} associated with an atomic number. If the atomic number does not exist, it returns a zero. | |
unsigned int | getTotalNumberOfSphericalFunctionsInCurrentProcessor () |
Returns the total number of total spherical functions indexed by {ilm} present in the current processor. If the atomic number does not exist, it returns a zero. | |
unsigned int | getMaximumNumberOfSphericalFunctions () |
Returns the maximum number of total spherical functions indexed by {ilm} across all atom Types present in atomNumbers vector. | |
void | getTotalAtomsAndNonLocalElementsInCurrentProcessor (unsigned int &totalAtomsInCurrentProcessor, unsigned int &totalNonLocalElements, std::vector< unsigned int > &numberCellsForEachAtom, std::vector< unsigned int > &numberCellsAccumNonLocalAtoms, std::vector< unsigned int > &iElemNonLocalToElemIndexMap) |
unsigned int | getTotalNumberOfRadialSphericalFunctions () |
Returns the total number of total radial-spherical functions indexed by {i} present in atomicNumbers list. | |
const std::map< std::pair< unsigned int, unsigned int >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > & | getSphericalFunctions () const |
Returns the shared_ptr of AtomCenteredSphericalFunctionBase associated with std::pair(atomic Number and lQuantumNo) | |
const std::vector< unsigned int > & | getAtomicNumbers () const |
Returns the vector of size Natoms of all atoms in system. | |
const std::vector< unsigned int > & | getAtomIdsInCurrentProcess () const |
Returns the atomIds of atoms present in current processor. | |
const unsigned int | getTotalSphericalFunctionIndexStart (unsigned int Znum, unsigned int alpha) |
Returns the startIndex of spherical Function alpha associated with atomic number Znum. | |
template<typename NumberType> | |
void | computeSparseStructure (std::shared_ptr< dftfe::basis::FEBasisOperations< NumberType, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const unsigned int quadratureIndex, const double cutOffVal=1.0E-8, const unsigned int cutOffType=0) |
void | setImageCoordinates (const std::vector< int > &imageIds, const std::vector< std::vector< double > > &periodicCoords) |
const std::vector< int > & | getAtomIdsInElement (unsigned int iElem) |
const std::map< unsigned int, std::vector< int > > & | getSparsityPattern () |
bool | atomSupportInElement (unsigned int iElem) |
void | getDataForSparseStructure (const std::map< unsigned int, std::vector< int > > &sparsityPattern, const std::vector< std::vector< dealii::CellId > > &elementIdsInAtomCompactSupport, const std::vector< std::vector< unsigned int > > &elementIndexesInAtomCompactSupport, const std::vector< unsigned int > &atomIdsInCurrentProcess, unsigned int numberElements) |
const unsigned int | getOffsetLocation (const unsigned int iAtom) |
Public Attributes | |
std::vector< std::vector< unsigned int > > | d_elementIndexesInAtomCompactSupport |
Private Attributes | |
std::vector< double > | d_atomCoords |
std::vector< unsigned int > | d_atomicNumbers |
std::map< unsigned int, std::vector< double > > | d_periodicImageCoord |
std::map< std::pair< unsigned int, unsigned int >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > | d_sphericalFunctionsContainer |
std::map< unsigned int, unsigned int > | d_numRadialSphericalFunctions |
std::map< unsigned int, unsigned int > | d_numSphericalFunctions |
std::map< unsigned int, std::vector< int > > | d_sparsityPattern |
std::vector< std::vector< dealii::CellId > > | d_elementIdsInAtomCompactSupport |
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > | d_elementOneFieldIteratorsInAtomCompactSupport |
std::vector< unsigned int > | d_AtomIdsInCurrentProcess |
std::vector< unsigned int > | d_offsetLocation |
std::vector< std::vector< int > > | d_AtomIdsInElement |
std::map< unsigned int, std::vector< unsigned int > > | d_totalSphericalFunctionIndexStart |
bool dftfe::AtomCenteredSphericalFunctionContainer::atomSupportInElement | ( | unsigned int | iElem | ) |
void dftfe::AtomCenteredSphericalFunctionContainer::computeSparseStructure | ( | std::shared_ptr< dftfe::basis::FEBasisOperations< NumberType, double, dftfe::utils::MemorySpace::HOST > > & | basisOperationsPtr, |
const unsigned int | quadratureIndex, | ||
const double | cutOffVal = 1.0E-8, | ||
const unsigned int | cutOffType = 0 ) |
const std::vector< double > & dftfe::AtomCenteredSphericalFunctionContainer::getAtomCoordinates | ( | ) | const |
Returns the cooridnates of atom present in domain.
const std::vector< unsigned int > & dftfe::AtomCenteredSphericalFunctionContainer::getAtomicNumbers | ( | ) | const |
Returns the vector of size Natoms of all atoms in system.
const std::vector< unsigned int > & dftfe::AtomCenteredSphericalFunctionContainer::getAtomIdsInCurrentProcess | ( | ) | const |
Returns the atomIds of atoms present in current processor.
const std::vector< int > & dftfe::AtomCenteredSphericalFunctionContainer::getAtomIdsInElement | ( | unsigned int | iElem | ) |
void dftfe::AtomCenteredSphericalFunctionContainer::getDataForSparseStructure | ( | const std::map< unsigned int, std::vector< int > > & | sparsityPattern, |
const std::vector< std::vector< dealii::CellId > > & | elementIdsInAtomCompactSupport, | ||
const std::vector< std::vector< unsigned int > > & | elementIndexesInAtomCompactSupport, | ||
const std::vector< unsigned int > & | atomIdsInCurrentProcess, | ||
unsigned int | numberElements ) |
unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getMaximumNumberOfSphericalFunctions | ( | ) |
Returns the maximum number of total spherical functions indexed by {ilm} across all atom Types present in atomNumbers vector.
unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getNumAtomCentersSize | ( | ) |
Returns the number of atoms present in domain.
const unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getOffsetLocation | ( | const unsigned int | iAtom | ) |
const std::map< unsigned int, std::vector< double > > & dftfe::AtomCenteredSphericalFunctionContainer::getPeriodicImageCoordinatesList | ( | ) | const |
Returns the map of atomId vs vector of image coordinates.
const std::map< unsigned int, std::vector< int > > & dftfe::AtomCenteredSphericalFunctionContainer::getSparsityPattern | ( | ) |
const std::map< std::pair< unsigned int, unsigned int >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > & dftfe::AtomCenteredSphericalFunctionContainer::getSphericalFunctions | ( | ) | const |
Returns the shared_ptr of AtomCenteredSphericalFunctionBase associated with std::pair(atomic Number and lQuantumNo)
void dftfe::AtomCenteredSphericalFunctionContainer::getTotalAtomsAndNonLocalElementsInCurrentProcessor | ( | unsigned int & | totalAtomsInCurrentProcessor, |
unsigned int & | totalNonLocalElements, | ||
std::vector< unsigned int > & | numberCellsForEachAtom, | ||
std::vector< unsigned int > & | numberCellsAccumNonLocalAtoms, | ||
std::vector< unsigned int > & | iElemNonLocalToElemIndexMap ) |
[out] | totalAtomsInCurrentProcessor | number of atoms in current processor based on compact support |
[out] | totalNonLocalElements | number of nonLocal elements in current processor |
[out] | numberCellsForEachAtom | number of cells associated which each atom in the current processor. vecot of size totalAtomsInCurrentProcessor |
[out] | numberCellsAccumNonLocalAtoms | number of cells accumulated till iatom in current processor. vector of size totalAtomsInCurrentProcessor |
unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getTotalNumberOfRadialSphericalFunctions | ( | ) |
Returns the total number of total radial-spherical functions indexed by {i} present in atomicNumbers list.
unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getTotalNumberOfRadialSphericalFunctionsPerAtom | ( | unsigned int | atomicNumber | ) |
Returns the he number of radial spherical functions indexed by {i} associated with an atomic number. If the atomic number does not exist, it returns a zero.
unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getTotalNumberOfSphericalFunctionsInCurrentProcessor | ( | ) |
Returns the total number of total spherical functions indexed by {ilm} present in the current processor. If the atomic number does not exist, it returns a zero.
unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getTotalNumberOfSphericalFunctionsPerAtom | ( | unsigned int | atomicNumber | ) |
Returns the he number of total spherical functions indexed by {ilm} associated with an atomic number. If the atomic number does not exist, it returns a zero.
const unsigned int dftfe::AtomCenteredSphericalFunctionContainer::getTotalSphericalFunctionIndexStart | ( | unsigned int | Znum, |
unsigned int | alpha ) |
Returns the startIndex of spherical Function alpha associated with atomic number Znum.
void dftfe::AtomCenteredSphericalFunctionContainer::init | ( | const std::vector< unsigned int > & | atomicNumbers, |
const std::map< std::pair< unsigned int, unsigned int >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > & | listOfSphericalFunctions ) |
Initialises the class with the atomicNumbers of various atoms and the AtomCenteredSphericalFn of various spherical functions. This function is only called once per run.
[in] | atomicNumbers | vector of size Natoms storing the Znumbers of various atoms present |
[in] | listOfSphericalFunctions | map of std::pain (Znum, l) to the sphericalFUnction class shared pointer. |
void dftfe::AtomCenteredSphericalFunctionContainer::initaliseCoordinates | ( | const std::vector< double > & | atomCoords, |
const std::vector< std::vector< double > > & | periodicCoords, | ||
const std::vector< int > & | imageIds ) |
Initialises the position of atoms, the image posisiton and image ids after every update of atom positions.
[in] | atomCoords | vector of size 3*Natoms storing the X,Y,Z coordiantes of atom in cell. |
[in] | periodicCoords | vector of vector storing the image coordinates |
[in] | imageIds | the image Id of image atoms present in periodicCoords input |
void dftfe::AtomCenteredSphericalFunctionContainer::setImageCoordinates | ( | const std::vector< int > & | imageIds, |
const std::vector< std::vector< double > > & | periodicCoords ) |
|
private |
|
private |
|
private |
|
private |
|
private |
std::vector<std::vector<unsigned int> > dftfe::AtomCenteredSphericalFunctionContainer::d_elementIndexesInAtomCompactSupport |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |