DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace > Class Template Reference

#include <AtomicCenteredNonLocalOperator.h>

Public Member Functions

 AtomicCenteredNonLocalOperator (std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr, std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, memorySpace > > basisOperatorPtr, std::shared_ptr< AtomCenteredSphericalFunctionContainer > atomCenteredSphericalFunctionContainer, const MPI_Comm &mpi_comm_parent, const bool memOptMode=false, const bool computeSphericalFnTimesX=true, const bool useGlobalCMatrix=false)
 
void initialiseOperatorActionOnX (unsigned int kPointIndex)
 Resizes various internal data members and selects the kpoint of interest.
 
void initialiseFlattenedDataStructure (unsigned int waveFunctionBlockSize, dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened)
 initialises the multivector object, waveFunctionBlockSize and resizes various internal data members.
 
void intitialisePartitionerKPointsAndComputeCMatrixEntries (const bool updateSparsity, const std::vector< double > &kPointWeights, const std::vector< double > &kPointCoordinates, std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperHostPtr, const unsigned int quadratureIndex)
 calls internal function: initialisePartitioner, initialiseKpoint and computeCMatrixEntries
 
template<typename ValueTypeSrc>
void copyPartitionerKPointsAndComputeCMatrixEntries (const bool updateSparsity, const std::vector< double > &kPointWeights, const std::vector< double > &kPointCoordinates, std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperHostPtr, const unsigned int quadratureIndex, const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueTypeSrc, memorySpace > > nonLocalOperatorSrc)
 calls internal function: initialisePartitioner, initialiseKpoint and computeCMatrixEntries
 
const std::vector< unsigned int > & getNonlocalElementToCellIdVector () const
 
unsigned int getTotalAtomInCurrentProcessor () const
 
const dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > & getFlattenedNonLocalCellDofIndexToProcessDofIndexMap () const
 
unsigned int getTotalNonLocalElementsInCurrentProcessor () const
 
unsigned int getTotalNonLocalEntriesCurrentProcessor () const
 
unsigned int getMaxSingleAtomEntries () const
 
bool atomSupportInElement (unsigned int iElem) const
 
unsigned int getGlobalDofAtomIdSphericalFnPair (const unsigned int atomId, const unsigned int alpha) const
 
unsigned int getLocalIdOfDistributedVec (const unsigned int globalId) const
 
std::vector< unsigned int > & getNonLocalElemIdToLocalElemIdMap () const
 
std::vector< unsigned int > & getAtomWiseNumberCellsInCompactSupport () const
 
std::vector< unsigned int > & getAtomWiseNumberCellsAccumulated () const
 
const std::vector< ValueType > & getAtomCenteredKpointIndexedSphericalFnQuadValues () const
 
const std::vector< ValueType > & getAtomCenteredKpointTimesSphericalFnTimesDistFromAtomQuadValues () const
 
const std::map< unsigned int, std::vector< unsigned int > > & getCellIdToAtomIdsLocalCompactSupportMap () const
 
const std::vector< unsigned int > & getNonTrivialSphericalFnsPerCell () const
 
const std::vector< unsigned int > & getNonTrivialSphericalFnsCellStartIndex () const
 
const unsigned int getTotalNonTrivialSphericalFnsOverAllCells () const
 
const std::vector< unsigned int > & getNonTrivialAllCellsSphericalFnAlphaToElemIdMap () const
 
const std::map< unsigned int, std::vector< unsigned int > > & getAtomIdToNonTrivialSphericalFnCellStartIndex () const
 Required in configurational forces. Cummulative sphercial Fn Id. The size is numCells in processor.
 
const std::vector< unsigned int > & getSphericalFnTimesVectorFlattenedVectorLocalIds () const
 Returns the Flattened vector of sphericalFunctionIDs in order of atomIDs of atoms in processor.
 
const std::vector< unsigned int > & getOwnedAtomIdsInCurrentProcessor () const
 
void computeCconjtransCMatrix (const unsigned int atomId, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperPtr, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &Dinverse, dftfe::utils::MemoryStorage< ValueType, dftfe::utils::MemorySpace::HOST > PconjtransposePmatrix)
 Computes C^{T}D^{-1}C at the global level for atomId. This is required in PAW.
 
void applyVOnCconjtransX (const CouplingStructure couplingtype, const dftfe::utils::MemoryStorage< ValueType, memorySpace > &couplingMatrix, dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened, const bool flagCopyResultsToMatrix=true, const unsigned int kPointIndex=0)
 compute the action of coupling matrix on sphericalFunctionKetTimesVectorParFlattened.
 
void copyBackFromDistributedVectorToLocalDataStructure (dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened, const dftfe::utils::MemoryStorage< double, memorySpace > &scalingVector)
 After AllReduce function is called this will copy to the nonLocalOperatorClassDatastructure.
 
void applyAllReduceOnCconjtransX (dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened, const bool skipComm=false)
 copies the results from internal member to sphericalFunctionKetTimesVectorParFlattened, on which ghost values are called. crucial operation for completion of the full CconjtranX on all cells
 
void applyCconjtransOnX (const ValueType *X, const std::pair< unsigned int, unsigned int > cellRange)
 computes the results of CconjtransX on the cells of interst specied by cellRange
 
void applyCconjtransOnX (const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &X)
 computes the results of CconjtransX on nodal X vector
 
const ValueType * getCconjtansXLocalDataStructure (const unsigned int iAtom) const
 Returns the pointer of CTX stored in HOST memory.
 
void applyVCconjtransOnX (const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &src, const unsigned int kPointIndex, const CouplingStructure couplingtype, const dftfe::utils::MemoryStorage< ValueType, memorySpace > &couplingMatrix, dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened, const bool flagScaleInternalMatrix=false)
 completes the VCconjX on nodal vector src. The src vector must have all ghost nodes and constraint nodes updated.
 
void applyCVCconjtransOnX (const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &src, const unsigned int kPointIndex, const CouplingStructure couplingtype, const dftfe::utils::MemoryStorage< ValueType, memorySpace > &couplingMatrix, dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened, dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &dst)
 completes the action of CVCconjtranspose on nodal vector src. The src vector must have all ghost nodes and contraint nodes updated.
 
void applyCOnVCconjtransX (ValueType *Xout, const std::pair< unsigned int, unsigned int > cellRange)
 adds the result of CVCtX onto Xout for both CPU and GPU calls
 
void applyCOnVCconjtransX (dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &Xout)
 adds the result of CVCtX onto Xout for both CPU and GPU calls
 
std::vector< ValueType > getCmatrixEntries (int kPointIndex, unsigned int atomId, int iElem) const
 
bool atomPresentInCellRange (const std::pair< unsigned int, unsigned int > cellRange) const
 
void paddingCouplingMatrix (const std::vector< ValueType > &entries, std::vector< ValueType > &entriesPadded, const CouplingStructure couplingtype)
 Called only for GPU runs where the coupling matrix has to be padded.
 
const std::vector< ValueType > & getCmatrixEntriesConjugate (const unsigned int chargeId, const unsigned int iElemComp) const
 Returns C matrix entries for chargeId and it compact support element Id.
 
const std::vector< ValueType > & getCmatrixEntriesTranspose (const unsigned int chargeId, const unsigned int iElemComp) const
 Returns C conj matrix entries for chargeId and it compact support element Id.
 
const std::vector< std::vector< dftfe::utils::MemoryStorage< ValueType, memorySpace > > > & getGlobalCMatrix () const
 Returns global C matrix of all atoms.
 

Protected Member Functions

void applyVCconjtransOnXCellLevel (const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &src, const unsigned int kPointIndex, const CouplingStructure couplingtype, const dftfe::utils::MemoryStorage< ValueType, memorySpace > &couplingMatrix, dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened, const bool flagScaleInternalMatrix=false)
 completes the VCconjX on nodal vector src. The src vector must have all ghost nodes and constraint nodes updated.
 
void applyVCconjtransOnXUsingGlobalC (const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &src, const unsigned int kPointIndex, const CouplingStructure couplingtype, const dftfe::utils::MemoryStorage< ValueType, memorySpace > &couplingMatrix, dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > &sphericalFunctionKetTimesVectorParFlattened, const bool flagScaleInternalMatrix=false)
 completes the VCconjX on nodal vector src using global C matrix. The global C matrix mush have been computed before. The src vector must have all ghost nodes and constraint nodes updated.
 

Protected Attributes

bool d_AllReduceCompleted
 
std::vector< double > d_kPointWeights
 
std::vector< double > d_kPointCoordinates
 
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
 
std::shared_ptr< AtomCenteredSphericalFunctionContainerd_atomCenteredSphericalFunctionContainer
 
std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > d_mpiPatternP2P
 
std::vector< unsigned int > d_numberCellsForEachAtom
 
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, memorySpace > > d_basisOperatorPtr
 
std::vector< ValueType > d_atomCenteredKpointIndexedSphericalFnQuadValues
 
std::vector< ValueType > d_atomCenteredKpointTimesSphericalFnTimesDistFromAtomQuadValues
 
std::map< unsigned int, std::vector< unsigned int > > d_cellIdToAtomIdsLocalCompactSupportMap
 map from cell number to set of non local atom ids (local numbering)
 
std::vector< unsigned int > d_nonTrivialSphericalFnPerCell
 vector of size num physical cells
 
std::vector< unsigned int > d_nonTrivialSphericalFnsCellStartIndex
 
std::vector< unsigned int > d_nonTrivialAllCellsSphericalFnAlphaToElemIdMap
 
std::map< unsigned int, std::vector< unsigned int > > d_atomIdToNonTrivialSphericalFnCellStartIndex
 map from local nonlocal atomid to vector over cells
 
unsigned int d_sumNonTrivialSphericalFnOverAllCells
 
std::vector< unsigned int > d_sphericalFnTimesVectorFlattenedVectorLocalIds
 
std::vector< distributedCPUVec< double > > d_SphericalFunctionKetTimesVectorPar
 
std::map< std::pair< unsigned int, unsigned int >, unsigned int > d_sphericalFunctionIdsNumberingMapCurrentProcess
 
std::vector< unsigned int > d_OwnedAtomIdsInCurrentProcessor
 
dealii::IndexSet d_locallyOwnedAtomCenteredFnIdsCurrentProcess
 
dealii::IndexSet d_ghostAtomCenteredFnIdsCurrentProcess
 
std::map< std::pair< unsigned int, unsigned int >, unsigned int > d_AtomCenteredFnIdsNumberingMapCurrentProcess
 
std::vector< std::vector< std::vector< dftfe::utils::MemoryStorage< ValueType, memorySpace > > > > d_CMatrixEntries
 
dealii::ConditionalOStream pcout
 
const MPI_Comm d_mpi_communicator
 
const unsigned int d_this_mpi_process
 
const unsigned int d_n_mpi_processes
 
dealii::IndexSet d_locallyOwnedSphericalFunctionIdsCurrentProcess
 
dealii::IndexSet d_ghostSphericalFunctionIdsCurrentProcess
 
unsigned int d_totalAtomsInCurrentProc
 
unsigned int d_totalNonlocalElems
 
unsigned int d_totalNonLocalEntries
 
unsigned int d_maxSingleAtomContribution
 
std::vector< unsigned int > d_numberCellsAccumNonLocalAtoms
 
dftfe::utils::MemoryStorage< unsigned int, memorySpace > d_iElemNonLocalToElemIndexMap
 
unsigned int d_numberNodesPerElement
 
unsigned int d_locallyOwnedCells
 
unsigned int d_numberWaveFunctions
 
unsigned int d_kPointIndex
 
bool d_memoryOptMode
 
bool d_isMallocCalled = false
 
std::vector< std::vector< std::vector< ValueType > > > d_CMatrixEntriesConjugate
 
std::vector< std::vector< std::vector< ValueType > > > d_CMatrixEntriesTranspose
 

Private Member Functions

void initKpoints (const std::vector< double > &kPointWeights, const std::vector< double > &kPointCoordinates)
 stores the d_kpointWeights, d_kpointCoordinates. Other data members regarding are computed from container data object
 
void initialisePartitioner ()
 creates the partitioner for the distributed vector based on sparsity patten from sphericalFn container.
 
void computeCMatrixEntries (std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr, const unsigned int quadratureIndex)
 computes the entries in C matrix for CPUs and GPUs. On GPUs the entries are copied to a flattened vector on device memory. Further on GPUs, various maps are created crucial for accessing and padding entries in Cmatrix flattened device.
 
template<typename ValueTypeSrc>
void copyCMatrixEntries (const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueTypeSrc, memorySpace > > nonLocalOperatorSrc, std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr, const unsigned int quadratureIndex)
 
template<typename ValueTypeSrc>
void copyGlobalCMatrix (const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueTypeSrc, memorySpace > > nonLocalOperatorSrc, std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr, const unsigned int quadratureIndex)
 
void computeGlobalCMatrixVector (std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperHostPtr)
 computes Global Cmatrix on HOST.
 

Private Attributes

std::map< unsigned int, dftfe::utils::MemoryStorage< ValueType, dftfe::utils::MemorySpace::HOST > > d_sphericalFnTimesWavefunMatrix
 
std::vector< dftfe::global_size_typed_flattenedNonLocalCellDofIndexToProcessDofIndexVector
 
dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > d_flattenedNonLocalCellDofIndexToProcessDofIndexMap
 
std::vector< unsigned int > d_nonlocalElemIdToCellIdVector
 
bool d_computeSphericalFnTimesX
 
bool d_useGlobalCMatrix
 
std::vector< unsigned int > d_atomStartIndexGlobal
 
unsigned int d_totalNumSphericalFunctionsGlobal
 
std::vector< std::vector< dftfe::utils::MemoryStorage< ValueType, memorySpace > > > d_CMatrixGlobal
 
std::set< unsigned int > d_setOfAtomicNumber
 
std::vector< unsigned int > d_mapAtomIdToSpeciesIndex
 
std::vector< unsigned int > d_mapiAtomToSpeciesIndex
 
std::vector< dftfe::utils::MemoryStorage< ValueType, memorySpace > > d_dotProductAtomicWaveInputWaveTemp
 
std::vector< unsigned int > d_mapIAtomicNumToDotProd
 
std::vector< unsigned int > d_mapiAtomToDotProd
 
unsigned int d_totalLocallyOwnedNodes
 
std::vector< unsigned int > d_mapiAtomTosphFuncWaveStart
 
std::map< unsigned int, std::vector< unsigned int > > d_listOfiAtomInSpecies
 

Constructor & Destructor Documentation

◆ AtomicCenteredNonLocalOperator()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::AtomicCenteredNonLocalOperator ( std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr,
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, memorySpace > > basisOperatorPtr,
std::shared_ptr< AtomCenteredSphericalFunctionContainer > atomCenteredSphericalFunctionContainer,
const MPI_Comm & mpi_comm_parent,
const bool memOptMode = false,
const bool computeSphericalFnTimesX = true,
const bool useGlobalCMatrix = false )

Member Function Documentation

◆ applyAllReduceOnCconjtransX()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyAllReduceOnCconjtransX ( dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened,
const bool skipComm = false )

copies the results from internal member to sphericalFunctionKetTimesVectorParFlattened, on which ghost values are called. crucial operation for completion of the full CconjtranX on all cells

Parameters
[in]sphericalFunctionKetTimesVectorParFlattenedmultivector to store results of CconjtransX which is initiliased using initialiseFlattenedVector call
[in]skip1flag for compute-communication overlap in ChFSI on GPUs
[in]skip2flag for compute-communication overlap in ChFSI on GPUs

◆ applyCconjtransOnX() [1/2]

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyCconjtransOnX ( const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & X)

computes the results of CconjtransX on nodal X vector

Parameters
[in]Xinput X nodal vector elements

◆ applyCconjtransOnX() [2/2]

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyCconjtransOnX ( const ValueType * X,
const std::pair< unsigned int, unsigned int > cellRange )

computes the results of CconjtransX on the cells of interst specied by cellRange

Parameters
[in]Xinput cell level vector
[in]cellRangestart and end element id in list of nonlocal elements

◆ applyCOnVCconjtransX() [1/2]

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyCOnVCconjtransX ( dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & Xout)

adds the result of CVCtX onto Xout for both CPU and GPU calls

Parameters
[out]XoutmemoryStorage object of size cells*numberOfNodex*BlockSize. Typical case holds the results of H_{loc}X
[in]cellRangestart and end element id in list of nonlocal elements

◆ applyCOnVCconjtransX() [2/2]

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyCOnVCconjtransX ( ValueType * Xout,
const std::pair< unsigned int, unsigned int > cellRange )

adds the result of CVCtX onto Xout for both CPU and GPU calls

Parameters
[out]XoutmemoryStorage object of size cells*numberOfNodex*BlockSize. Typical case holds the results of H_{loc}X
[in]cellRangestart and end element id in list of nonlocal elements

◆ applyCVCconjtransOnX()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyCVCconjtransOnX ( const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & src,
const unsigned int kPointIndex,
const CouplingStructure couplingtype,
const dftfe::utils::MemoryStorage< ValueType, memorySpace > & couplingMatrix,
dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened,
dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & dst )

completes the action of CVCconjtranspose on nodal vector src. The src vector must have all ghost nodes and contraint nodes updated.

Parameters
[in]srcinput nodal vector on which operator acts on.
[in]kPointIndexkPoint of interst for current operation
[in]couplingtypestructure of coupling matrix
[in]couplingMatrixentires of the coupling matrix V in CVCconjtrans
[in]sphericalFunctionKetTimesVectorParFlattenedmultivector to store results of CconjtransX which is initiliased using initialiseFlattenedVector call
[out]dstoutput nodal vector where the results of the operator is copied into.

◆ applyVCconjtransOnX()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyVCconjtransOnX ( const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & src,
const unsigned int kPointIndex,
const CouplingStructure couplingtype,
const dftfe::utils::MemoryStorage< ValueType, memorySpace > & couplingMatrix,
dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened,
const bool flagScaleInternalMatrix = false )

completes the VCconjX on nodal vector src. The src vector must have all ghost nodes and constraint nodes updated.

Parameters
[in]srcinput nodal vector on which operator acts on.
[in]kPointIndexkPoint of interest for current operation
[in]couplingtypestructure of coupling matrix
[in]couplingMatrixentries of the coupling matrix V in CVCconjtrans. Ensure the coupling matrix is padded
[out]sphericalFunctionKetTimesVectorParFlattenedmultivector to store results of CconjtransX which is initiliased using initialiseFlattenedVector call

◆ applyVCconjtransOnXCellLevel()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyVCconjtransOnXCellLevel ( const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & src,
const unsigned int kPointIndex,
const CouplingStructure couplingtype,
const dftfe::utils::MemoryStorage< ValueType, memorySpace > & couplingMatrix,
dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened,
const bool flagScaleInternalMatrix = false )
protected

completes the VCconjX on nodal vector src. The src vector must have all ghost nodes and constraint nodes updated.

Parameters
[in]srcinput nodal vector on which operator acts on.
[in]kPointIndexkPoint of interest for current operation
[in]couplingtypestructure of coupling matrix
[in]couplingMatrixentries of the coupling matrix V in CVCconjtrans. Ensure the coupling matrix is padded
[out]sphericalFunctionKetTimesVectorParFlattenedmultivector to store results of CconjtransX which is initiliased using initialiseFlattenedVector call

◆ applyVCconjtransOnXUsingGlobalC()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyVCconjtransOnXUsingGlobalC ( const dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & src,
const unsigned int kPointIndex,
const CouplingStructure couplingtype,
const dftfe::utils::MemoryStorage< ValueType, memorySpace > & couplingMatrix,
dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened,
const bool flagScaleInternalMatrix = false )
protected

completes the VCconjX on nodal vector src using global C matrix. The global C matrix mush have been computed before. The src vector must have all ghost nodes and constraint nodes updated.

Parameters
[in]srcinput nodal vector on which operator acts on.
[in]kPointIndexkPoint of interest for current operation
[in]couplingtypestructure of coupling matrix
[in]couplingMatrixentries of the coupling matrix V in CVCconjtrans. Ensure the coupling matrix is padded
[out]sphericalFunctionKetTimesVectorParFlattenedmultivector to store results of CconjtransX which is initiliased using initialiseFlattenedVector call

◆ applyVOnCconjtransX()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::applyVOnCconjtransX ( const CouplingStructure couplingtype,
const dftfe::utils::MemoryStorage< ValueType, memorySpace > & couplingMatrix,
dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened,
const bool flagCopyResultsToMatrix = true,
const unsigned int kPointIndex = 0 )

compute the action of coupling matrix on sphericalFunctionKetTimesVectorParFlattened.

Parameters
[in]couplingtypestructure of coupling matrix
[in]couplingMatrixentires of the coupling matrix V in CVCconjtrans. Ensure that the coupling matrix is padded. Refer to ONCVclass for template
[out]sphericalFunctionKetTimesVectorParFlattenedmultivector to store results of CconjtransX which is initiliased using initialiseFlattenedVector call. The results are stored in sphericalFunctionKetTimesVectorParFlattened or internal data member based on flagCopyResultsToMatrix.
[in]flagCopyResultsToMatrixflag to confirm whether to scal the multivector sphericalFunctionKetTimesVectorParFlattened or store results in internal data member.

◆ atomPresentInCellRange()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
bool dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::atomPresentInCellRange ( const std::pair< unsigned int, unsigned int > cellRange) const

◆ atomSupportInElement()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
bool dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::atomSupportInElement ( unsigned int iElem) const

◆ computeCconjtransCMatrix()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::computeCconjtransCMatrix ( const unsigned int atomId,
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperPtr,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & Dinverse,
dftfe::utils::MemoryStorage< ValueType, dftfe::utils::MemorySpace::HOST > PconjtransposePmatrix )

Computes C^{T}D^{-1}C at the global level for atomId. This is required in PAW.

◆ computeCMatrixEntries()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::computeCMatrixEntries ( std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr,
const unsigned int quadratureIndex )
private

computes the entries in C matrix for CPUs and GPUs. On GPUs the entries are copied to a flattened vector on device memory. Further on GPUs, various maps are created crucial for accessing and padding entries in Cmatrix flattened device.

Parameters
[in]basisOperationsPtrHOST FEBasisOperations shared_ptr required to indetify the element ids and quad points
[in]quadratureIndexquadrature index for sampling the spherical function. Quadrature Index is used to reinit basisOperationsPtr

◆ computeGlobalCMatrixVector()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::computeGlobalCMatrixVector ( std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr,
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperHostPtr )
private

computes Global Cmatrix on HOST.

Parameters
[in]basisOperationsPtrHOST FEBasisOperations shared_ptr required to indetify the element ids and quad points
[in]BLASWrapperHostPtrHOST BLASWrapper

◆ copyBackFromDistributedVectorToLocalDataStructure()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::copyBackFromDistributedVectorToLocalDataStructure ( dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened,
const dftfe::utils::MemoryStorage< double, memorySpace > & scalingVector )

After AllReduce function is called this will copy to the nonLocalOperatorClassDatastructure.

◆ copyCMatrixEntries()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
template<typename ValueTypeSrc>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::copyCMatrixEntries ( const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueTypeSrc, memorySpace > > nonLocalOperatorSrc,
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr,
const unsigned int quadratureIndex )
private

◆ copyGlobalCMatrix()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
template<typename ValueTypeSrc>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::copyGlobalCMatrix ( const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueTypeSrc, memorySpace > > nonLocalOperatorSrc,
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr,
const unsigned int quadratureIndex )
private

◆ copyPartitionerKPointsAndComputeCMatrixEntries()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
template<typename ValueTypeSrc>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::copyPartitionerKPointsAndComputeCMatrixEntries ( const bool updateSparsity,
const std::vector< double > & kPointWeights,
const std::vector< double > & kPointCoordinates,
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr,
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperHostPtr,
const unsigned int quadratureIndex,
const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueTypeSrc, memorySpace > > nonLocalOperatorSrc )

calls internal function: initialisePartitioner, initialiseKpoint and computeCMatrixEntries

Parameters
[in]updateSparsityflag on whether the sparstiy patten was updated, hence the partitioner is updated.
[in]kPointWeightsstd::vector<double> of size number of kPoints
[out]kPointCoordinatesstd::vector<double> of kPoint coordinates
[in]basisOperationsPtrHOST FEBasisOperations shared_ptr required to indetify the element ids and quad points
[in]BLASWrapperHostPtrCPU blasWrapperPtr, used for xcopy calls
[in]quadratureIndexquadrature index for sampling the spherical function. Quadrature Index is used to reinit basisOperationsPtr
[in]nonLocalOperatorSrcThe source nonLocalOpertor from where the CMatrix and partitioner is copied. Generally, it is of higher precision.

◆ getAtomCenteredKpointIndexedSphericalFnQuadValues()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< ValueType > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getAtomCenteredKpointIndexedSphericalFnQuadValues ( ) const

◆ getAtomCenteredKpointTimesSphericalFnTimesDistFromAtomQuadValues()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< ValueType > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getAtomCenteredKpointTimesSphericalFnTimesDistFromAtomQuadValues ( ) const

◆ getAtomIdToNonTrivialSphericalFnCellStartIndex()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::map< unsigned int, std::vector< unsigned int > > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getAtomIdToNonTrivialSphericalFnCellStartIndex ( ) const

Required in configurational forces. Cummulative sphercial Fn Id. The size is numCells in processor.

◆ getAtomWiseNumberCellsAccumulated()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getAtomWiseNumberCellsAccumulated ( ) const

◆ getAtomWiseNumberCellsInCompactSupport()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getAtomWiseNumberCellsInCompactSupport ( ) const

◆ getCconjtansXLocalDataStructure()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const ValueType * dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getCconjtansXLocalDataStructure ( const unsigned int iAtom) const

Returns the pointer of CTX stored in HOST memory.

Parameters
[in]iAtomatomIndex in the list of atoms with support in the current processor. NOTE!! One must be careful here

◆ getCellIdToAtomIdsLocalCompactSupportMap()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::map< unsigned int, std::vector< unsigned int > > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getCellIdToAtomIdsLocalCompactSupportMap ( ) const

◆ getCmatrixEntries()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector< ValueType > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getCmatrixEntries ( int kPointIndex,
unsigned int atomId,
int iElem ) const

◆ getCmatrixEntriesConjugate()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< ValueType > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getCmatrixEntriesConjugate ( const unsigned int chargeId,
const unsigned int iElemComp ) const

Returns C matrix entries for chargeId and it compact support element Id.

◆ getCmatrixEntriesTranspose()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< ValueType > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getCmatrixEntriesTranspose ( const unsigned int chargeId,
const unsigned int iElemComp ) const

Returns C conj matrix entries for chargeId and it compact support element Id.

◆ getFlattenedNonLocalCellDofIndexToProcessDofIndexMap()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getFlattenedNonLocalCellDofIndexToProcessDofIndexMap ( ) const

◆ getGlobalCMatrix()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< std::vector< dftfe::utils::MemoryStorage< ValueType, memorySpace > > > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getGlobalCMatrix ( ) const

Returns global C matrix of all atoms.

◆ getGlobalDofAtomIdSphericalFnPair()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getGlobalDofAtomIdSphericalFnPair ( const unsigned int atomId,
const unsigned int alpha ) const

◆ getLocalIdOfDistributedVec()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getLocalIdOfDistributedVec ( const unsigned int globalId) const

◆ getMaxSingleAtomEntries()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getMaxSingleAtomEntries ( ) const

◆ getNonlocalElementToCellIdVector()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getNonlocalElementToCellIdVector ( ) const

◆ getNonLocalElemIdToLocalElemIdMap()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getNonLocalElemIdToLocalElemIdMap ( ) const

◆ getNonTrivialAllCellsSphericalFnAlphaToElemIdMap()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getNonTrivialAllCellsSphericalFnAlphaToElemIdMap ( ) const

◆ getNonTrivialSphericalFnsCellStartIndex()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getNonTrivialSphericalFnsCellStartIndex ( ) const

◆ getNonTrivialSphericalFnsPerCell()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getNonTrivialSphericalFnsPerCell ( ) const

◆ getOwnedAtomIdsInCurrentProcessor()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getOwnedAtomIdsInCurrentProcessor ( ) const

◆ getSphericalFnTimesVectorFlattenedVectorLocalIds()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const std::vector< unsigned int > & dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getSphericalFnTimesVectorFlattenedVectorLocalIds ( ) const

Returns the Flattened vector of sphericalFunctionIDs in order of atomIDs of atoms in processor.

◆ getTotalAtomInCurrentProcessor()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getTotalAtomInCurrentProcessor ( ) const

◆ getTotalNonLocalElementsInCurrentProcessor()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getTotalNonLocalElementsInCurrentProcessor ( ) const

◆ getTotalNonLocalEntriesCurrentProcessor()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getTotalNonLocalEntriesCurrentProcessor ( ) const

◆ getTotalNonTrivialSphericalFnsOverAllCells()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::getTotalNonTrivialSphericalFnsOverAllCells ( ) const

◆ initialiseFlattenedDataStructure()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::initialiseFlattenedDataStructure ( unsigned int waveFunctionBlockSize,
dftfe::linearAlgebra::MultiVector< ValueType, memorySpace > & sphericalFunctionKetTimesVectorParFlattened )

initialises the multivector object, waveFunctionBlockSize and resizes various internal data members.

Parameters
[in]waveFunctionBlockSizesets the wavefunction block size for the action of the nonlocal operator.
[out]sphericalFunctionKetTimesVectorParFlattened,themultivector that is initialised based on blocksize and partitioner.

◆ initialiseOperatorActionOnX()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::initialiseOperatorActionOnX ( unsigned int kPointIndex)

Resizes various internal data members and selects the kpoint of interest.

Parameters
[in]kPointIndexspecifies the k-point of interest

◆ initialisePartitioner()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::initialisePartitioner ( )
private

creates the partitioner for the distributed vector based on sparsity patten from sphericalFn container.

Parameters
[in]basisOperationsPtrHOST FEBasisOperations shared_ptr required to indetify the element ids and quad points.

◆ initKpoints()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::initKpoints ( const std::vector< double > & kPointWeights,
const std::vector< double > & kPointCoordinates )
private

stores the d_kpointWeights, d_kpointCoordinates. Other data members regarding are computed from container data object

Parameters
[in]kPointWeightsstd::vector<double> of size number of kPoints
[out]kPointCoordinatesstd::vector<double> of kPoint coordinates

◆ intitialisePartitionerKPointsAndComputeCMatrixEntries()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::intitialisePartitionerKPointsAndComputeCMatrixEntries ( const bool updateSparsity,
const std::vector< double > & kPointWeights,
const std::vector< double > & kPointCoordinates,
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > basisOperationsPtr,
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperHostPtr,
const unsigned int quadratureIndex )

calls internal function: initialisePartitioner, initialiseKpoint and computeCMatrixEntries

Parameters
[in]updateSparsityflag on whether the sparstiy patten was updated, hence the partitioner is updated.
[in]kPointWeightsstd::vector<double> of size number of kPoints
[out]kPointCoordinatesstd::vector<double> of kPoint coordinates
[in]basisOperationsPtrHOST FEBasisOperations shared_ptr required to indetify the element ids and quad points
[in]quadratureIndexquadrature index for sampling the spherical function. Quadrature Index is used to reinit basisOperationsPtr

◆ paddingCouplingMatrix()

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
void dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::paddingCouplingMatrix ( const std::vector< ValueType > & entries,
std::vector< ValueType > & entriesPadded,
const CouplingStructure couplingtype )

Called only for GPU runs where the coupling matrix has to be padded.

Parameters
[in]entriesCOupling matrix entries without padding in the atomId order
[out]entriesPaddedPadding of coupling matrix entries
[in]couplingtypeDetermines the dimension of entriesPadded and the padding mechanism elements

Member Data Documentation

◆ d_AllReduceCompleted

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
bool dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_AllReduceCompleted
protected

◆ d_AtomCenteredFnIdsNumberingMapCurrentProcess

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::map<std::pair<unsigned int, unsigned int>, unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_AtomCenteredFnIdsNumberingMapCurrentProcess
protected

◆ d_atomCenteredKpointIndexedSphericalFnQuadValues

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<ValueType> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_atomCenteredKpointIndexedSphericalFnQuadValues
protected

◆ d_atomCenteredKpointTimesSphericalFnTimesDistFromAtomQuadValues

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<ValueType> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_atomCenteredKpointTimesSphericalFnTimesDistFromAtomQuadValues
protected

◆ d_atomCenteredSphericalFunctionContainer

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<AtomCenteredSphericalFunctionContainer> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_atomCenteredSphericalFunctionContainer
protected

◆ d_atomIdToNonTrivialSphericalFnCellStartIndex

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::map<unsigned int, std::vector<unsigned int> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_atomIdToNonTrivialSphericalFnCellStartIndex
protected

map from local nonlocal atomid to vector over cells

◆ d_atomStartIndexGlobal

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_atomStartIndexGlobal
private

◆ d_basisOperatorPtr

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::shared_ptr< dftfe::basis::FEBasisOperations<dataTypes::number, double, memorySpace> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_basisOperatorPtr
protected

◆ d_BLASWrapperPtr

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_BLASWrapperPtr
protected

◆ d_cellIdToAtomIdsLocalCompactSupportMap

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::map<unsigned int, std::vector<unsigned int> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_cellIdToAtomIdsLocalCompactSupportMap
protected

map from cell number to set of non local atom ids (local numbering)

◆ d_CMatrixEntries

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<std::vector< std::vector<dftfe::utils::MemoryStorage<ValueType, memorySpace> > > > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_CMatrixEntries
protected

◆ d_CMatrixEntriesConjugate

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<std::vector<std::vector<ValueType> > > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_CMatrixEntriesConjugate
protected

◆ d_CMatrixEntriesTranspose

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<std::vector<std::vector<ValueType> > > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_CMatrixEntriesTranspose
protected

◆ d_CMatrixGlobal

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector< std::vector<dftfe::utils::MemoryStorage<ValueType, memorySpace> > > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_CMatrixGlobal
private

◆ d_computeSphericalFnTimesX

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
bool dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_computeSphericalFnTimesX
private

◆ d_dotProductAtomicWaveInputWaveTemp

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<dftfe::utils::MemoryStorage<ValueType, memorySpace> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_dotProductAtomicWaveInputWaveTemp
private

◆ d_flattenedNonLocalCellDofIndexToProcessDofIndexMap

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<dftfe::global_size_type, memorySpace> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_flattenedNonLocalCellDofIndexToProcessDofIndexMap
private

◆ d_flattenedNonLocalCellDofIndexToProcessDofIndexVector

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<dftfe::global_size_type> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_flattenedNonLocalCellDofIndexToProcessDofIndexVector
private

◆ d_ghostAtomCenteredFnIdsCurrentProcess

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_ghostAtomCenteredFnIdsCurrentProcess
protected

◆ d_ghostSphericalFunctionIdsCurrentProcess

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_ghostSphericalFunctionIdsCurrentProcess
protected

◆ d_iElemNonLocalToElemIndexMap

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<unsigned int, memorySpace> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_iElemNonLocalToElemIndexMap
protected

◆ d_isMallocCalled

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
bool dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_isMallocCalled = false
protected

◆ d_kPointCoordinates

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_kPointCoordinates
protected

◆ d_kPointIndex

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_kPointIndex
protected

◆ d_kPointWeights

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_kPointWeights
protected

◆ d_listOfiAtomInSpecies

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::map<unsigned int, std::vector<unsigned int> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_listOfiAtomInSpecies
private

◆ d_locallyOwnedAtomCenteredFnIdsCurrentProcess

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_locallyOwnedAtomCenteredFnIdsCurrentProcess
protected

◆ d_locallyOwnedCells

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_locallyOwnedCells
protected

◆ d_locallyOwnedSphericalFunctionIdsCurrentProcess

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_locallyOwnedSphericalFunctionIdsCurrentProcess
protected

◆ d_mapAtomIdToSpeciesIndex

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_mapAtomIdToSpeciesIndex
private

◆ d_mapIAtomicNumToDotProd

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_mapIAtomicNumToDotProd
private

◆ d_mapiAtomToDotProd

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_mapiAtomToDotProd
private

◆ d_mapiAtomToSpeciesIndex

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_mapiAtomToSpeciesIndex
private

◆ d_mapiAtomTosphFuncWaveStart

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_mapiAtomTosphFuncWaveStart
private

◆ d_maxSingleAtomContribution

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_maxSingleAtomContribution
protected

◆ d_memoryOptMode

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
bool dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_memoryOptMode
protected

◆ d_mpi_communicator

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const MPI_Comm dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_mpi_communicator
protected

◆ d_mpiPatternP2P

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::shared_ptr< const utils::mpi::MPIPatternP2P<dftfe::utils::MemorySpace::HOST> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_mpiPatternP2P
protected

◆ d_n_mpi_processes

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_n_mpi_processes
protected

◆ d_nonlocalElemIdToCellIdVector

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_nonlocalElemIdToCellIdVector
private

◆ d_nonTrivialAllCellsSphericalFnAlphaToElemIdMap

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_nonTrivialAllCellsSphericalFnAlphaToElemIdMap
protected

◆ d_nonTrivialSphericalFnPerCell

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_nonTrivialSphericalFnPerCell
protected

vector of size num physical cells

◆ d_nonTrivialSphericalFnsCellStartIndex

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_nonTrivialSphericalFnsCellStartIndex
protected

vector of size num physical cell with starting index for each cell for the above array

◆ d_numberCellsAccumNonLocalAtoms

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_numberCellsAccumNonLocalAtoms
protected

◆ d_numberCellsForEachAtom

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_numberCellsForEachAtom
protected

◆ d_numberNodesPerElement

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_numberNodesPerElement
protected

◆ d_numberWaveFunctions

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_numberWaveFunctions
protected

◆ d_OwnedAtomIdsInCurrentProcessor

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_OwnedAtomIdsInCurrentProcessor
protected

◆ d_setOfAtomicNumber

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::set<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_setOfAtomicNumber
private

◆ d_sphericalFnTimesVectorFlattenedVectorLocalIds

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_sphericalFnTimesVectorFlattenedVectorLocalIds
protected

◆ d_sphericalFnTimesWavefunMatrix

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::map< unsigned int, dftfe::utils::MemoryStorage<ValueType, dftfe::utils::MemorySpace::HOST> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_sphericalFnTimesWavefunMatrix
private

◆ d_sphericalFunctionIdsNumberingMapCurrentProcess

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::map<std::pair<unsigned int, unsigned int>, unsigned int> dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_sphericalFunctionIdsNumberingMapCurrentProcess
protected

◆ d_SphericalFunctionKetTimesVectorPar

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
std::vector<distributedCPUVec<double> > dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_SphericalFunctionKetTimesVectorPar
protected

◆ d_sumNonTrivialSphericalFnOverAllCells

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_sumNonTrivialSphericalFnOverAllCells
protected

◆ d_this_mpi_process

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_this_mpi_process
protected

◆ d_totalAtomsInCurrentProc

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_totalAtomsInCurrentProc
protected

◆ d_totalLocallyOwnedNodes

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_totalLocallyOwnedNodes
private

◆ d_totalNonlocalElems

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_totalNonlocalElems
protected

◆ d_totalNonLocalEntries

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_totalNonLocalEntries
protected

◆ d_totalNumSphericalFunctionsGlobal

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_totalNumSphericalFunctionsGlobal
private

◆ d_useGlobalCMatrix

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
bool dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::d_useGlobalCMatrix
private

◆ pcout

template<typename ValueType, dftfe::utils::MemorySpace memorySpace>
dealii::ConditionalOStream dftfe::AtomicCenteredNonLocalOperator< ValueType, memorySpace >::pcout
protected

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