20#ifndef DFTFE_ATOMCENTEREDSPHERICALFUNCTIONCONTAINERBASE_H
21#define DFTFE_ATOMCENTEREDSPHERICALFUNCTIONCONTAINERBASE_H
47 init(
const std::vector<dftfe::uInt> &atomicNumbers,
48 const std::map<std::pair<dftfe::uInt, dftfe::uInt>,
49 std::shared_ptr<AtomCenteredSphericalFunctionBase>>
50 &listOfSphericalFunctions);
61 const std::vector<std::vector<double>> &periodicCoords,
62 const std::vector<dftfe::Int> &imageIds);
75 const std::vector<double> &
81 const std::map<dftfe::uInt, std::vector<double>> &
126 std::vector<dftfe::uInt> &numberCellsForEachAtom,
127 std::vector<dftfe::uInt> &numberCellsAccumNonLocalAtoms,
128 std::vector<dftfe::uInt> &iElemNonLocalToElemIndexMap);
139 const std::map<std::pair<dftfe::uInt, dftfe::uInt>,
140 std::shared_ptr<AtomCenteredSphericalFunctionBase>> &
145 const std::vector<dftfe::uInt> &
150 const std::vector<dftfe::uInt> &
161 template <
typename NumberType>
170 const double cutOffVal = 1.0E-8,
177 const std::vector<std::vector<double>> &periodicCoords);
181 const std::vector<dftfe::Int> &
184 const std::map<dftfe::uInt, std::vector<dftfe::Int>> &
192 const std::map<
dftfe::uInt, std::vector<dftfe::Int>> &sparsityPattern,
193 const std::vector<std::vector<dealii::CellId>>
194 &elementIdsInAtomCompactSupport,
195 const std::vector<std::vector<dftfe::uInt>>
196 &elementIndexesInAtomCompactSupport,
197 const std::vector<dftfe::uInt> &atomIdsInCurrentProcess,
223 std::map<std::pair<dftfe::uInt, dftfe::uInt>,
224 std::shared_ptr<AtomCenteredSphericalFunctionBase>>
238 std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator>>
243 std::map<dftfe::uInt, std::vector<dftfe::uInt>>
Definition AtomCenteredSphericalFunctionContainer.h:37
void init(const std::vector< dftfe::uInt > &atomicNumbers, const std::map< std::pair< dftfe::uInt, dftfe::uInt >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > &listOfSphericalFunctions)
Initialises the class with the atomicNumbers of various atoms and the AtomCenteredSphericalFn of vari...
bool atomSupportInElement(dftfe::uInt iElem)
void computeSparseStructure(std::shared_ptr< dftfe::basis::FEBasisOperations< NumberType, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt quadratureIndex, const double cutOffVal=1.0E-8, const dftfe::uInt cutOffType=0)
std::map< dftfe::uInt, dftfe::uInt > d_numRadialSphericalFunctions
Definition AtomCenteredSphericalFunctionContainer.h:227
std::vector< std::vector< dealii::CellId > > d_elementIdsInAtomCompactSupport
Definition AtomCenteredSphericalFunctionContainer.h:235
const std::vector< dftfe::Int > & getAtomIdsInElement(dftfe::uInt iElem)
const std::vector< dftfe::uInt > & getAtomIdsInCurrentProcess() const
Returns the atomIds of atoms present in current processor.
void setImageCoordinates(const std::vector< dftfe::Int > &imageIds, const std::vector< std::vector< double > > &periodicCoords)
void initaliseCoordinates(const std::vector< double > &atomCoords, const std::vector< std::vector< double > > &periodicCoords, const std::vector< dftfe::Int > &imageIds)
Initialises the position of atoms, the image posisiton and image ids after every update of atom posit...
std::vector< dftfe::uInt > d_atomicNumbers
Definition AtomCenteredSphericalFunctionContainer.h:213
std::vector< dftfe::uInt > d_offsetLocation
Definition AtomCenteredSphericalFunctionContainer.h:241
const std::map< dftfe::uInt, std::vector< dftfe::Int > > & getSparsityPattern()
void getTotalAtomsAndNonLocalElementsInCurrentProcessor(dftfe::uInt &totalAtomsInCurrentProcessor, dftfe::uInt &totalNonLocalElements, std::vector< dftfe::uInt > &numberCellsForEachAtom, std::vector< dftfe::uInt > &numberCellsAccumNonLocalAtoms, std::vector< dftfe::uInt > &iElemNonLocalToElemIndexMap)
const std::vector< dftfe::uInt > & getAtomicNumbers() const
Returns the vector of size Natoms of all atoms in system.
dftfe::uInt getTotalNumberOfSphericalFunctionsPerAtom(dftfe::uInt atomicNumber)
Returns the he number of total spherical functions indexed by {ilm} associated with an atomic number....
std::map< dftfe::uInt, std::vector< dftfe::Int > > d_sparsityPattern
Definition AtomCenteredSphericalFunctionContainer.h:233
std::vector< dftfe::uInt > d_AtomIdsInCurrentProcess
Definition AtomCenteredSphericalFunctionContainer.h:240
dftfe::uInt getTotalNumberOfSphericalFunctionsInCurrentProcessor()
Returns the total number of total spherical functions indexed by {ilm} present in the current process...
void getDataForSparseStructure(const std::map< dftfe::uInt, std::vector< dftfe::Int > > &sparsityPattern, const std::vector< std::vector< dealii::CellId > > &elementIdsInAtomCompactSupport, const std::vector< std::vector< dftfe::uInt > > &elementIndexesInAtomCompactSupport, const std::vector< dftfe::uInt > &atomIdsInCurrentProcess, dftfe::uInt numberElements)
const dftfe::uInt getOffsetLocation(const dftfe::uInt iAtom)
std::map< dftfe::uInt, std::vector< dftfe::uInt > > d_totalSphericalFunctionIndexStart
Definition AtomCenteredSphericalFunctionContainer.h:244
const std::map< std::pair< dftfe::uInt, dftfe::uInt >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > & getSphericalFunctions() const
Returns the shared_ptr of AtomCenteredSphericalFunctionBase associated with std::pair(atomic Number a...
const dftfe::uInt getTotalSphericalFunctionIndexStart(dftfe::uInt Znum, dftfe::uInt alpha)
Returns the startIndex of spherical Function alpha associated with atomic number Znum.
std::vector< double > d_atomCoords
Definition AtomCenteredSphericalFunctionContainer.h:208
const std::vector< double > & getAtomCoordinates() const
Returns the cooridnates of atom present in domain.
dftfe::uInt getTotalNumberOfRadialSphericalFunctions()
Returns the total number of total radial-spherical functions indexed by {i} present in atomicNumbers ...
std::map< std::pair< dftfe::uInt, dftfe::uInt >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_sphericalFunctionsContainer
Definition AtomCenteredSphericalFunctionContainer.h:225
dftfe::uInt getTotalNumberOfRadialSphericalFunctionsPerAtom(dftfe::uInt atomicNumber)
Returns the he number of radial spherical functions indexed by {i} associated with an atomic number....
std::map< dftfe::uInt, std::vector< double > > d_periodicImageCoord
Definition AtomCenteredSphericalFunctionContainer.h:218
std::map< dftfe::uInt, dftfe::uInt > d_numSphericalFunctions
Definition AtomCenteredSphericalFunctionContainer.h:230
std::vector< std::vector< dftfe::uInt > > d_elementIndexesInAtomCompactSupport
Definition AtomCenteredSphericalFunctionContainer.h:174
dftfe::uInt getMaximumNumberOfSphericalFunctions()
Returns the maximum number of total spherical functions indexed by {ilm} across all atom Types presen...
dftfe::uInt getNumAtomCentersSize()
Returns the number of atoms present in domain.
std::vector< std::vector< dftfe::Int > > d_AtomIdsInElement
Definition AtomCenteredSphericalFunctionContainer.h:242
const std::map< dftfe::uInt, std::vector< double > > & getPeriodicImageCoordinatesList() const
Returns the map of atomId vs vector of image coordinates.
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_elementOneFieldIteratorsInAtomCompactSupport
Definition AtomCenteredSphericalFunctionContainer.h:239
Definition FEBasisOperations.h:85
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10