DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
|
#include <FEBasisOperations.h>
Public Member Functions | |
FEBasisOperations (std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr) | |
Constructor. | |
~FEBasisOperations ()=default | |
Default Destructor. | |
void | clear () |
Clears the FEBasisOperations internal storage. | |
void | init (dealii::MatrixFree< 3, ValueTypeBasisData > &matrixFreeData, std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector, const dftfe::uInt &dofHandlerID, const std::vector< dftfe::uInt > &quadratureID, const std::vector< UpdateFlags > updateFlags) |
fills required data structures for the given dofHandlerID | |
template<dftfe::utils::MemorySpace memorySpaceSrc> | |
void | init (const FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpaceSrc > &basisOperationsSrc) |
fills required data structures from another FEBasisOperations object | |
void | reinit (const dftfe::uInt &vecBlockSize, const dftfe::uInt &cellBlockSize, const dftfe::uInt &quadratureID, const bool isResizeTempStorageForInerpolation=true, const bool isResizeTempStorageForCellMatrices=false) |
sets internal variables and optionally resizes internal temp storage for interpolation operations | |
dftfe::utils::MemoryStorage< dftfe::uInt, dftfe::utils::MemorySpace::HOST > & | getFlattenedMapsHost () |
void | initializeIndexMaps () |
Initializes indexset maps from process level indices to cell level indices for a single vector, also initializes cell index to cellid map. | |
void | initializeFlattenedIndexMaps () |
Initializes indexset maps from process level indices to cell level indices for multivectors. | |
void | initializeConstraints () |
Initializes the constraintMatrixInfo object. | |
void | reinitializeConstraints (std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector) |
Reinitializes the constraintMatrixInfo object. | |
void | initializeMPIPattern () |
Constructs the MPIPatternP2P object. | |
void | initializeShapeFunctionAndJacobianData () |
Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype. | |
void | initializeShapeFunctionAndJacobianBasisData () |
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype. | |
void | computeCellStiffnessMatrix (const dftfe::uInt quadratureID, const dftfe::uInt cellsBlockSize, const bool basisType=false, const bool ceoffType=true) |
Computes the cell-level stiffness matrix. | |
void | computeCellMassMatrix (const dftfe::uInt quadratureID, const dftfe::uInt cellsBlockSize, const bool basisType=false, const bool ceoffType=true) |
void | computeWeightedCellMassMatrix (const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellMassMatrix) const |
void | computeWeightedCellNjGradNiMatrix (const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiMatrix) const |
void | computeWeightedCellNjGradNiPlusNiGradNjMatrix (const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const |
void | computeWeightedCellNjGradNiMinusNiGradNjMatrix (const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const |
void | computeWeightedCellStiffnessMatrix (const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellStiffnessMatrix) const |
void | computeInverseSqrtMassVector (const bool basisType=true, const bool ceoffType=false) |
void | computeStiffnessVector (const bool basisType=true, const bool ceoffType=false) |
Computes the stiffness matrix \grad Ni . \grad Ni. | |
void | resizeTempStorage (const bool isResizeTempStorageForInerpolation, const bool isResizeTempStorageForCellMatrices) |
Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in reinit. | |
dftfe::uInt | nQuadsPerCell () const |
Number of quadrature points per cell for the quadratureID set in reinit. | |
dftfe::uInt | nVectors () const |
Number of vectors set in reinit. | |
dftfe::uInt | nDofsPerCell () const |
Number of DoFs per cell for the dofHandlerID set in init. | |
dftfe::uInt | nCells () const |
Number of locally owned cells on the current processor. | |
dftfe::uInt | nRelaventDofs () const |
Number of DoFs on the current processor, locally owned + ghosts. | |
dftfe::uInt | nOwnedDofs () const |
Number of locally owned DoFs on the current processor. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & | shapeFunctionData (bool transpose=false) const |
Shape function values at quadrature points. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & | shapeFunctionGradientData (bool transpose=false) const |
Shape function gradient values at quadrature points. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & | inverseJacobians () const |
Inverse Jacobian matrices, for cartesian cells returns the diagonal elements of the inverse Jacobian matrices for each cell, for affine cells returns the 3x3 inverse Jacobians for each cell otherwise returns the 3x3 inverse Jacobians at each quad point for each cell. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & | JxW () const |
determinant of Jacobian times the quadrature weight at each quad point for each cell. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > & | quadPoints () const |
quad point coordinates for each cell. | |
const auto & | shapeFunctionBasisData (bool transpose=false) const |
Shape function values at quadrature points in ValueTypeBasisData. | |
const auto & | shapeFunctionGradientBasisData (bool transpose=false) const |
Shape function gradient values at quadrature points in ValueTypeBasisData. | |
const auto & | inverseJacobiansBasisData () const |
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of the inverse Jacobian matrices for each cell, for affine cells returns the 3x3 inverse Jacobians for each cell otherwise returns the 3x3 inverse Jacobians at each quad point for each cell. | |
const auto & | JxWBasisData () const |
determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each quad point for each cell. | |
const auto & | cellStiffnessMatrix () const |
Cell level stiffness matrix in ValueTypeBasisCoeff. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | cellStiffnessMatrixBasisData () const |
Cell level stiffness matrix in ValueTypeBasisData. | |
const auto & | cellMassMatrix () const |
Cell level mass matrix in ValueTypeBasisCoeff. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | cellMassMatrixBasisData () const |
Cell level mass matrix in ValueTypeBasisData. | |
const auto & | cellInverseSqrtMassVector () const |
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff. | |
const auto & | cellInverseMassVector () const |
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff. | |
const auto & | cellSqrtMassVector () const |
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff. | |
const auto & | cellMassVector () const |
Cell level diagonal mass matrix in ValueTypeBasisCoeff. | |
const auto & | inverseSqrtMassVector () const |
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff. | |
const auto & | inverseMassVector () const |
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff. | |
const auto & | sqrtMassVector () const |
sqrt diagonal mass matrix in ValueTypeBasisCoeff | |
const auto & | massVector () const |
diagonal mass matrix in ValueTypeBasisCoeff | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | cellInverseSqrtMassVectorBasisData () const |
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | cellInverseMassVectorBasisData () const |
Cell level inverse diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & | cellInverseMassVectorBasisDataSinglePrec () const |
Cell level inverse diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | cellSqrtMassVectorBasisData () const |
Cell level sqrt diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | cellMassVectorBasisData () const |
Cell level diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | inverseSqrtMassVectorBasisData () const |
Inverse sqrt diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | inverseMassVectorBasisData () const |
Inverse diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & | inverseMassVectorBasisDataSinglePrec () const |
Inverse diagonal mass matrix in ValueTypeBasisData. | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | sqrtMassVectorBasisData () const |
sqrt diagonal mass matrix in ValueTypeBasisData | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | massVectorBasisData () const |
diagonal mass matrix in ValueTypeBasisData | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | stiffnessVectorBasisData () const |
diagonal stiffness matrix in ValueTypeBasisData | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | inverseStiffnessVectorBasisData () const |
diagonal inverse stiffness matrix in ValueTypeBasisData | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | inverseSqrtStiffnessVectorBasisData () const |
diagonal inverse sqrt stiffness matrix in ValueTypeBasisData | |
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | sqrtStiffnessVectorBasisData () const |
diagonal sqrt stiffness matrix in ValueTypeBasisData | |
dftfe::uInt | cellsTypeFlag () const |
returns 2 if all cells on current processor are Cartesian, 1 if all cells on current processor are affine and 0 otherwise. | |
dealii::CellId | cellID (const dftfe::uInt iElem) const |
returns the deal.ii cellID corresponing to given cell Index. | |
dealii::DoFHandler< 3 >::active_cell_iterator | getCellIterator (const dftfe::uInt iElem) const |
returns the deal.ii cell_iterator corresponing to given cell Index. | |
dftfe::uInt | cellIndex (const dealii::CellId cellid) const |
returns the cell index corresponding to given deal.ii cellID. | |
void | createMultiVector (const dftfe::uInt blocksize, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector) const |
Creates a multivector. | |
void | createMultiVectorSinglePrec (const dftfe::uInt blocksize, dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > &multiVector) const |
Creates a multivector. | |
void | createScratchMultiVectors (const dftfe::uInt vecBlockSize, const dftfe::uInt numMultiVecs=1) const |
Creates scratch multivectors. | |
void | createScratchMultiVectorsSinglePrec (const dftfe::uInt vecBlockSize, const dftfe::uInt numMultiVecs=1) const |
Creates single precision scratch multivectors. | |
void | clearScratchMultiVectors () const |
Clears scratch multivectors. | |
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | getMultiVector (const dftfe::uInt vecBlockSize, const dftfe::uInt index=0) const |
Gets scratch multivectors. | |
dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > & | getMultiVectorSinglePrec (const dftfe::uInt vecBlockSize, const dftfe::uInt index=0) const |
Gets single precision scratch multivectors. | |
void | distribute (dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector, dftfe::uInt constraintIndex=std::numeric_limits< dftfe::uInt >::max()) const |
Apply constraints on given multivector. | |
const dealii::MatrixFree< 3, ValueTypeBasisData > & | matrixFreeData () const |
Return the underlying deal.II matrixfree object. | |
const dealii::DoFHandler< 3 > & | getDofHandler () const |
Return the underlying deal.II dofhandler object. | |
void | interpolate (dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients=NULL) const |
Interpolate process level nodal data to cell level quadrature data. | |
void | integrateWithBasis (ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > &mapQuadIdToProcId) const |
Integrate cell level quadrature data times shape functions to process level nodal data. | |
void | extractToCellNodalData (dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr) const |
Get cell level nodal data from process level nodal data. | |
void | accumulateFromCellNodalData (const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData) const |
Accumulate cell level nodal data into process level nodal data. | |
void | interpolateKernel (const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const |
Interpolate process level nodal data to cell level quadrature data. | |
void | interpolateKernel (const ValueTypeBasisCoeff *nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const |
Interpolate cell level nodal data to cell level quadrature data. | |
void | integrateWithBasisKernel (const ValueTypeBasisCoeff *quadratureValues, const ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > &mapQuadIdToProcId, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const |
Integrate cell level quadrature data times shape functions to process level nodal data. | |
void | extractToCellNodalDataKernel (const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const |
Get cell level nodal data from process level nodal data. | |
void | accumulateFromCellNodalDataKernel (const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const |
Accumulate cell level nodal data into process level nodal data. | |
Protected Attributes | |
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > | tempCellNodalData |
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > | tempQuadratureGradientsData |
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > | tempQuadratureGradientsDataNonAffine |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > | tempCellGradientsBlock |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > | tempCellGradientsBlock2 |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > | tempCellValuesBlock |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > | tempCellMatrixBlock |
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > | zeroIndexVec |
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > | tempCellValuesBlockCoeff |
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > | d_BLASWrapperPtr |
dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::FEBasisOperations | ( | std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > | BLASWrapperPtr | ) |
Constructor.
|
default |
Default Destructor.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::accumulateFromCellNodalData | ( | const ValueTypeBasisCoeff * | cellNodalDataPtr, |
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData ) const |
Accumulate cell level nodal data into process level nodal data.
[in] | cellNodalDataPtr | Cell level nodal values, indexed by [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec]. |
[out] | nodalData | process level nodal data. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::accumulateFromCellNodalDataKernel | ( | const ValueTypeBasisCoeff * | cellNodalDataPtr, |
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData, | ||
const std::pair< dftfe::uInt, dftfe::uInt > | cellRange ) const |
Accumulate cell level nodal data into process level nodal data.
[in] | cellNodalDataPtr | Cell level nodal values, indexed by [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec]. |
[out] | nodalData | process level nodal data. |
[in] | cellRange | the range of cells for which extraction has to be done. |
dealii::CellId dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellID | ( | const dftfe::uInt | iElem | ) | const |
returns the deal.ii cellID corresponing to given cell Index.
[in] | iElem | cell Index |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellIndex | ( | const dealii::CellId | cellid | ) | const |
returns the cell index corresponding to given deal.ii cellID.
[in] | iElem | cell Index |
|
inline |
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellInverseMassVectorBasisData | ( | ) | const |
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellInverseMassVectorBasisDataSinglePrec | ( | ) | const |
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
|
inline |
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellInverseSqrtMassVectorBasisData | ( | ) | const |
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData.
|
inline |
Cell level mass matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellMassMatrixBasisData | ( | ) | const |
Cell level mass matrix in ValueTypeBasisData.
|
inline |
Cell level diagonal mass matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellMassVectorBasisData | ( | ) | const |
Cell level diagonal mass matrix in ValueTypeBasisData.
|
inline |
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellSqrtMassVectorBasisData | ( | ) | const |
Cell level sqrt diagonal mass matrix in ValueTypeBasisData.
|
inline |
Cell level stiffness matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellStiffnessMatrixBasisData | ( | ) | const |
Cell level stiffness matrix in ValueTypeBasisData.
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::cellsTypeFlag | ( | ) | const |
returns 2 if all cells on current processor are Cartesian, 1 if all cells on current processor are affine and 0 otherwise.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::clear | ( | ) |
Clears the FEBasisOperations internal storage.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::clearScratchMultiVectors | ( | ) | const |
Clears scratch multivectors.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeCellMassMatrix | ( | const dftfe::uInt | quadratureID, |
const dftfe::uInt | cellsBlockSize, | ||
const bool | basisType = false, | ||
const bool | ceoffType = true ) |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeCellStiffnessMatrix | ( | const dftfe::uInt | quadratureID, |
const dftfe::uInt | cellsBlockSize, | ||
const bool | basisType = false, | ||
const bool | ceoffType = true ) |
Computes the cell-level stiffness matrix.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeInverseSqrtMassVector | ( | const bool | basisType = true, |
const bool | ceoffType = false ) |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeStiffnessVector | ( | const bool | basisType = true, |
const bool | ceoffType = false ) |
Computes the stiffness matrix \grad Ni . \grad Ni.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeWeightedCellMassMatrix | ( | const std::pair< dftfe::uInt, dftfe::uInt > | cellRangeTotal, |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weights, | ||
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weightedCellMassMatrix ) const |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeWeightedCellNjGradNiMatrix | ( | const std::pair< dftfe::uInt, dftfe::uInt > | cellRangeTotal, |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weights, | ||
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weightedCellNjGradNiMatrix ) const |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeWeightedCellNjGradNiMinusNiGradNjMatrix | ( | const std::pair< dftfe::uInt, dftfe::uInt > | cellRangeTotal, |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weights, | ||
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weightedCellNjGradNiPlusNiGradNjMatrix ) const |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeWeightedCellNjGradNiPlusNiGradNjMatrix | ( | const std::pair< dftfe::uInt, dftfe::uInt > | cellRangeTotal, |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weights, | ||
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weightedCellNjGradNiPlusNiGradNjMatrix ) const |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::computeWeightedCellStiffnessMatrix | ( | const std::pair< dftfe::uInt, dftfe::uInt > | cellRangeTotal, |
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weights, | ||
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & | weightedCellStiffnessMatrix ) const |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::createMultiVector | ( | const dftfe::uInt | blocksize, |
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | multiVector ) const |
Creates a multivector.
[in] | blocksize | Number of vectors in the multivector. |
[out] | multiVector | the created multivector. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::createMultiVectorSinglePrec | ( | const dftfe::uInt | blocksize, |
dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > & | multiVector ) const |
Creates a multivector.
[in] | blocksize | Number of vectors in the multivector. |
[out] | multiVector | the created multivector. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::createScratchMultiVectors | ( | const dftfe::uInt | vecBlockSize, |
const dftfe::uInt | numMultiVecs = 1 ) const |
Creates scratch multivectors.
[in] | vecBlockSize | Number of vectors in the multivector. |
[out] | numMultiVecs | number of scratch multivectors needed with this vecBlockSize. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::createScratchMultiVectorsSinglePrec | ( | const dftfe::uInt | vecBlockSize, |
const dftfe::uInt | numMultiVecs = 1 ) const |
Creates single precision scratch multivectors.
[in] | vecBlockSize | Number of vectors in the multivector. |
[out] | numMultiVecs | number of scratch multivectors needed with this vecBlockSize. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::distribute | ( | dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | multiVector, |
dftfe::uInt | constraintIndex = std::numeric_limits< dftfe::uInt >::max() ) const |
Apply constraints on given multivector.
[in,out] | multiVector | the given multivector. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::extractToCellNodalData | ( | dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData, |
ValueTypeBasisCoeff * | cellNodalDataPtr ) const |
Get cell level nodal data from process level nodal data.
[in] | nodalData | process level nodal data, the multivector should already have ghost data and constraints should have been applied. |
[out] | cellNodalDataPtr | Cell level nodal values, indexed by [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec]. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::extractToCellNodalDataKernel | ( | const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData, |
ValueTypeBasisCoeff * | cellNodalDataPtr, | ||
const std::pair< dftfe::uInt, dftfe::uInt > | cellRange ) const |
Get cell level nodal data from process level nodal data.
[in] | nodalData | process level nodal data, the multivector should already have ghost data and constraints should have been applied. |
[out] | cellNodalDataPtr | Cell level nodal values, indexed by [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec]. |
[in] | cellRange | the range of cells for which extraction has to be done. |
dealii::DoFHandler< 3 >::active_cell_iterator dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::getCellIterator | ( | const dftfe::uInt | iElem | ) | const |
returns the deal.ii cell_iterator corresponing to given cell Index.
[in] | iElem | cell Index |
const dealii::DoFHandler< 3 > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::getDofHandler | ( | ) | const |
Return the underlying deal.II dofhandler object.
dftfe::utils::MemoryStorage< dftfe::uInt, dftfe::utils::MemorySpace::HOST > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::getFlattenedMapsHost | ( | ) |
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::getMultiVector | ( | const dftfe::uInt | vecBlockSize, |
const dftfe::uInt | index = 0 ) const |
Gets scratch multivectors.
[in] | vecBlockSize | Number of vectors in the multivector. |
[out] | numMultiVecs | index of the multivector among those with the same vecBlockSize. |
dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::getMultiVectorSinglePrec | ( | const dftfe::uInt | vecBlockSize, |
const dftfe::uInt | index = 0 ) const |
Gets single precision scratch multivectors.
[in] | vecBlockSize | Number of vectors in the multivector. |
[out] | numMultiVecs | index of the multivector among those with the same vecBlockSize. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::init | ( | const FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpaceSrc > & | basisOperationsSrc | ) |
fills required data structures from another FEBasisOperations object
[in] | basisOperationsSrc | Source FEBasisOperations object. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::init | ( | dealii::MatrixFree< 3, ValueTypeBasisData > & | matrixFreeData, |
std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > & | constraintsVector, | ||
const dftfe::uInt & | dofHandlerID, | ||
const std::vector< dftfe::uInt > & | quadratureID, | ||
const std::vector< UpdateFlags > | updateFlags ) |
fills required data structures for the given dofHandlerID
[in] | matrixFreeData | MatrixFree object. |
[in] | constraintsVector | std::vector of AffineConstraints, should be the same vector which was passed for the construction of the given MatrixFree object. |
[in] | dofHandlerID | dofHandler index to be used for getting data from the MatrixFree object. |
[in] | quadratureID | std::vector of quadratureIDs to be used, should be the same IDs which were used during the construction of the given MatrixFree object. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::initializeConstraints | ( | ) |
Initializes the constraintMatrixInfo object.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::initializeFlattenedIndexMaps | ( | ) |
Initializes indexset maps from process level indices to cell level indices for multivectors.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::initializeIndexMaps | ( | ) |
Initializes indexset maps from process level indices to cell level indices for a single vector, also initializes cell index to cellid map.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::initializeMPIPattern | ( | ) |
Constructs the MPIPatternP2P object.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::initializeShapeFunctionAndJacobianBasisData | ( | ) |
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::initializeShapeFunctionAndJacobianData | ( | ) |
Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::integrateWithBasis | ( | ValueTypeBasisCoeff * | quadratureValues, |
ValueTypeBasisCoeff * | quadratureGradients, | ||
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData, | ||
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > & | mapQuadIdToProcId ) const |
Integrate cell level quadrature data times shape functions to process level nodal data.
[in] | quadratureValues | Cell level quadrature values, indexed by [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[in] | quadratureGradients | Cell level quadrature gradients, indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[out] | nodalData | process level nodal data. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::integrateWithBasisKernel | ( | const ValueTypeBasisCoeff * | quadratureValues, |
const ValueTypeBasisCoeff * | quadratureGradients, | ||
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData, | ||
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > & | mapQuadIdToProcId, | ||
const std::pair< dftfe::uInt, dftfe::uInt > | cellRange ) const |
Integrate cell level quadrature data times shape functions to process level nodal data.
[in] | quadratureValues | Cell level quadrature values, indexed by [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[in] | quadratureGradients | Cell level quadrature gradients, indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[out] | nodalData | process level nodal data. |
[in] | cellRange | the range of cells for which integration has to be done. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::interpolate | ( | dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData, |
ValueTypeBasisCoeff * | quadratureValues, | ||
ValueTypeBasisCoeff * | quadratureGradients = NULL ) const |
Interpolate process level nodal data to cell level quadrature data.
[in] | nodalData | process level nodal data, the multivector should already have ghost data and constraints should have been applied. |
[out] | quadratureValues | Cell level quadrature values, indexed by [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[out] | quadratureGradients | Cell level quadrature gradients, indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::interpolateKernel | ( | const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & | nodalData, |
ValueTypeBasisCoeff * | quadratureValues, | ||
ValueTypeBasisCoeff * | quadratureGradients, | ||
const std::pair< dftfe::uInt, dftfe::uInt > | cellRange ) const |
Interpolate process level nodal data to cell level quadrature data.
[in] | nodalData | process level nodal data, the multivector should already have ghost data and constraints should have been applied. |
[out] | quadratureValues | Cell level quadrature values, indexed by [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[out] | quadratureGradients | Cell level quadrature gradients, indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[in] | cellRange | the range of cells for which interpolation has to be done. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::interpolateKernel | ( | const ValueTypeBasisCoeff * | nodalData, |
ValueTypeBasisCoeff * | quadratureValues, | ||
ValueTypeBasisCoeff * | quadratureGradients, | ||
const std::pair< dftfe::uInt, dftfe::uInt > | cellRange ) const |
Interpolate cell level nodal data to cell level quadrature data.
[in] | nodalData | cell level nodal data, the multivector should already have ghost data and constraints should have been applied. |
[out] | quadratureValues | Cell level quadrature values, indexed by [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[out] | quadratureGradients | Cell level quadrature gradients, indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec]. |
[in] | cellRange | the range of cells for which interpolation has to be done. |
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::inverseJacobians | ( | ) | const |
Inverse Jacobian matrices, for cartesian cells returns the diagonal elements of the inverse Jacobian matrices for each cell, for affine cells returns the 3x3 inverse Jacobians for each cell otherwise returns the 3x3 inverse Jacobians at each quad point for each cell.
|
inline |
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of the inverse Jacobian matrices for each cell, for affine cells returns the 3x3 inverse Jacobians for each cell otherwise returns the 3x3 inverse Jacobians at each quad point for each cell.
|
inline |
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::inverseMassVectorBasisData | ( | ) | const |
Inverse diagonal mass matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::inverseMassVectorBasisDataSinglePrec | ( | ) | const |
Inverse diagonal mass matrix in ValueTypeBasisData.
|
inline |
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::inverseSqrtMassVectorBasisData | ( | ) | const |
Inverse sqrt diagonal mass matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::inverseSqrtStiffnessVectorBasisData | ( | ) | const |
diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::inverseStiffnessVectorBasisData | ( | ) | const |
diagonal inverse stiffness matrix in ValueTypeBasisData
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::JxW | ( | ) | const |
determinant of Jacobian times the quadrature weight at each quad point for each cell.
|
inline |
determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each quad point for each cell.
|
inline |
diagonal mass matrix in ValueTypeBasisCoeff
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::massVectorBasisData | ( | ) | const |
diagonal mass matrix in ValueTypeBasisData
const dealii::MatrixFree< 3, ValueTypeBasisData > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::matrixFreeData | ( | ) | const |
Return the underlying deal.II matrixfree object.
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::nCells | ( | ) | const |
Number of locally owned cells on the current processor.
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::nDofsPerCell | ( | ) | const |
Number of DoFs per cell for the dofHandlerID set in init.
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::nOwnedDofs | ( | ) | const |
Number of locally owned DoFs on the current processor.
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::nQuadsPerCell | ( | ) | const |
Number of quadrature points per cell for the quadratureID set in reinit.
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::nRelaventDofs | ( | ) | const |
Number of DoFs on the current processor, locally owned + ghosts.
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::nVectors | ( | ) | const |
Number of vectors set in reinit.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::quadPoints | ( | ) | const |
quad point coordinates for each cell.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::reinit | ( | const dftfe::uInt & | vecBlockSize, |
const dftfe::uInt & | cellBlockSize, | ||
const dftfe::uInt & | quadratureID, | ||
const bool | isResizeTempStorageForInerpolation = true, | ||
const bool | isResizeTempStorageForCellMatrices = false ) |
sets internal variables and optionally resizes internal temp storage for interpolation operations
[in] | vecBlockSize | block size to used for operations on vectors, this has to be set to the exact value before any such operations are called. |
[in] | cellBlockSize | block size to used for cells, this has to be set to a value greater than or equal to the required value before any such operations are called |
[in] | quadratureID | Quadrature index to be used. |
[in] | isResizeTempStorage | whether to resize internal tempstorage. |
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::reinitializeConstraints | ( | std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > & | constraintsVector | ) |
Reinitializes the constraintMatrixInfo object.
void dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::resizeTempStorage | ( | const bool | isResizeTempStorageForInerpolation, |
const bool | isResizeTempStorageForCellMatrices ) |
Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in reinit.
|
inline |
Shape function values at quadrature points in ValueTypeBasisData.
[in] | transpose | if false the the data is indexed as [iQuad * d_nDofsPerCell + iNode] and if true it is indexed as [iNode * d_nQuadsPerCell + iQuad]. |
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::shapeFunctionData | ( | bool | transpose = false | ) | const |
Shape function values at quadrature points.
[in] | transpose | if false the the data is indexed as [iQuad * d_nDofsPerCell + iNode] and if true it is indexed as [iNode * d_nQuadsPerCell + iQuad]. |
|
inline |
Shape function gradient values at quadrature points in ValueTypeBasisData.
[in] | transpose | if false the the data is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell + iNode * d_nQuadsPerCell + iQuad]. |
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::shapeFunctionGradientData | ( | bool | transpose = false | ) | const |
Shape function gradient values at quadrature points.
[in] | transpose | if false the the data is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell + iNode * d_nQuadsPerCell + iQuad]. |
|
inline |
sqrt diagonal mass matrix in ValueTypeBasisCoeff
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::sqrtMassVectorBasisData | ( | ) | const |
sqrt diagonal mass matrix in ValueTypeBasisData
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::sqrtStiffnessVectorBasisData | ( | ) | const |
diagonal sqrt stiffness matrix in ValueTypeBasisData
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::stiffnessVectorBasisData | ( | ) | const |
diagonal stiffness matrix in ValueTypeBasisData
bool dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::areAllCellsAffine |
bool dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::areAllCellsCartesian |
|
protected |
dftfe::utils::MemoryStorage<dftfe::uInt, dftfe::utils::MemorySpace::HOST> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellDofIndexToProcessDofIndexMap |
std::map<dealii::CellId, dftfe::uInt> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellIdToCellIndexMap |
std::vector<dealii::CellId> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellIndexToCellIdMap |
std::vector<dealii::DoFHandler<3>::active_cell_iterator> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellIndexToCellIteratorMap |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseMassVectorBasisType |
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType<ValueTypeBasisData>::type, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseMassVectorBasisTypeSinglePrec |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseMassVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseSqrtMassVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseSqrtMassVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseSqrtStiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseSqrtStiffnessVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseStiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellInverseStiffnessVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellMassMatrixBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellMassMatrixCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellMassVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellMassVectorCoeffType |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellsBlockSize |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellSqrtMassVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellSqrtMassVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellSqrtStiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellSqrtStiffnessVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellStiffnessMatrixBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellStiffnessMatrixCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellStiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_cellStiffnessVectorCoeffType |
std::vector<dftUtils::constraintMatrixInfo<memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_constraintInfo |
std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>* dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_constraintsVector |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_dofHandlerID |
dftfe::utils::MemoryStorage<dftfe::uInt, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_flattenedCellDofIndexToProcessDofIndexMap |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseJacobianBasisData |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseJacobianData |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseMassVectorBasisType |
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType<ValueTypeBasisData>::type, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseMassVectorBasisTypeSinglePrec |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseMassVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseSqrtMassVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseSqrtMassVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseSqrtStiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseSqrtStiffnessVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseStiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_inverseStiffnessVectorCoeffType |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_JxWBasisData |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_JxWData |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_locallyOwnedSize |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_localSize |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_massVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_massVectorCoeffType |
const dealii::MatrixFree<3, ValueTypeBasisData>* dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_matrixFreeDataPtr |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_nCells |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_nDofsPerCell |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_nOMPThreads |
std::vector<dftfe::uInt> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_nQuadsPerCell |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_nVectors |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisData, dftfe::utils::MemorySpace::HOST> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_quadPoints |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_quadratureID |
std::vector<dftfe::uInt> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_quadratureIDsVector |
dftfe::uInt dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_quadratureIndex |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionBasisData |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionBasisDataTranspose |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionData |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionDataTranspose |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionGradientBasisData |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionGradientBasisDataTranspose |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionGradientData |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionGradientDataInternalLayout |
std::map<dftfe::uInt, dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_shapeFunctionGradientDataTranspose |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_sqrtMassVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_sqrtMassVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_sqrtStiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_sqrtStiffnessVectorCoeffType |
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_stiffnessVectorBasisType |
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_stiffnessVectorCoeffType |
std::vector<UpdateFlags> dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::d_updateFlags |
std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace> > dftfe::basis::FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpace >::mpiPatternP2P |
|
mutable |
|
mutable |
|
mutableprotected |
|
protected |
|
protected |
|
mutableprotected |
|
protected |
|
mutableprotected |
|
protected |
|
protected |
|
mutableprotected |