18#ifndef dftfeFEBasisOperations_h
19#define dftfeFEBasisOperations_h
51 return static_cast<UpdateFlags>(
static_cast<unsigned int>(f1) |
52 static_cast<unsigned int>(f2));
67 return static_cast<UpdateFlags>(
static_cast<unsigned int>(f1) &
68 static_cast<unsigned int>(f2));
80 template <
typename ValueTypeBasisCoeff,
81 typename ValueTypeBasisData,
96 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
133 std::vector<
const dealii::AffineConstraints<ValueTypeBasisData> *>
135 const unsigned int & dofHandlerID,
136 const std::vector<unsigned int> &quadratureID,
137 const std::vector<UpdateFlags> updateFlags);
143 template <dftfe::utils::MemorySpace memorySpaceSrc>
147 memorySpaceSrc> &basisOperationsSrc);
162 const unsigned int &cellBlockSize,
163 const unsigned int &quadratureID,
164 const bool isResizeTempStorageForInerpolation =
true,
165 const bool isResizeTempStorageForCellMatrices =
false);
197 std::vector<
const dealii::AffineConstraints<ValueTypeBasisData> *>
224 const unsigned int cellsBlockSize,
225 const bool basisType =
false,
226 const bool ceoffType =
true);
230 const unsigned int cellsBlockSize,
231 const bool basisType =
false,
232 const bool ceoffType =
true);
236 const std::pair<unsigned int, unsigned int> cellRangeTotal,
239 &weightedCellMassMatrix)
const;
243 const std::pair<unsigned int, unsigned int> cellRangeTotal,
246 &weightedCellNjGradNiMatrix)
const;
250 const std::pair<unsigned int, unsigned int> cellRangeTotal,
253 &weightedCellNjGradNiPlusNiGradNjMatrix)
const;
257 const bool ceoffType =
false);
265 const bool ceoffType =
false);
272 const bool isResizeTempStorageForCellMatrices);
360 if constexpr (std::is_same<ValueTypeBasisCoeff,
361 ValueTypeBasisData>::value)
385 if constexpr (std::is_same<ValueTypeBasisCoeff,
386 ValueTypeBasisData>::value)
413 if constexpr (std::is_same<ValueTypeBasisCoeff,
414 ValueTypeBasisData>::value)
435 if constexpr (std::is_same<ValueTypeBasisCoeff,
436 ValueTypeBasisData>::value)
452 if constexpr (std::is_same<ValueTypeBasisCoeff,
453 ValueTypeBasisData>::value)
477 if constexpr (std::is_same<ValueTypeBasisCoeff,
478 ValueTypeBasisData>::value)
502 if constexpr (std::is_same<ValueTypeBasisCoeff,
503 ValueTypeBasisData>::value)
519 if constexpr (std::is_same<ValueTypeBasisCoeff,
520 ValueTypeBasisData>::value)
536 if constexpr (std::is_same<ValueTypeBasisCoeff,
537 ValueTypeBasisData>::value)
553 if constexpr (std::is_same<ValueTypeBasisCoeff,
554 ValueTypeBasisData>::value)
571 if constexpr (std::is_same<ValueTypeBasisCoeff,
572 ValueTypeBasisData>::value)
590 if constexpr (std::is_same<ValueTypeBasisCoeff,
591 ValueTypeBasisData>::value)
608 if constexpr (std::is_same<ValueTypeBasisCoeff,
609 ValueTypeBasisData>::value)
625 if constexpr (std::is_same<ValueTypeBasisCoeff,
626 ValueTypeBasisData>::value)
745 dealii::DoFHandler<3>::active_cell_iterator
762 const unsigned int blocksize,
773 const unsigned int blocksize,
776 memorySpace> &multiVector)
const;
786 const unsigned int numMultiVecs = 1)
const;
796 const unsigned int vecBlockSize,
797 const unsigned int numMultiVecs = 1)
const;
813 const unsigned int index = 0)
const;
825 const unsigned int index = 0)
const;
833 memorySpace> &multiVector,
834 unsigned int constraintIndex =
835 std::numeric_limits<unsigned int>::max())
const;
842 const dealii::MatrixFree<3, ValueTypeBasisData> &
848 const dealii::DoFHandler<3> &
855 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
861 std::map<
unsigned int,
868 std::vector<dealii::DoFHandler<3>::active_cell_iterator>
871 std::map<
unsigned int,
874 std::map<
unsigned int,
877 std::map<
unsigned int,
880 std::map<
unsigned int,
883 std::map<
unsigned int,
886 std::map<
unsigned int,
889 std::map<
unsigned int,
893 std::map<
unsigned int,
896 std::map<
unsigned int,
899 std::map<
unsigned int,
902 std::map<
unsigned int,
905 std::map<
unsigned int,
908 std::map<
unsigned int,
1024 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
1040 memorySpace> &nodalData,
1041 ValueTypeBasisCoeff *quadratureValues,
1042 ValueTypeBasisCoeff *quadratureGradients = NULL)
const;
1056 ValueTypeBasisCoeff *quadratureValues,
1057 ValueTypeBasisCoeff *quadratureGradients,
1061 &mapQuadIdToProcId)
const;
1074 ValueTypeBasisCoeff *cellNodalDataPtr)
const;
1084 const ValueTypeBasisCoeff *cellNodalDataPtr,
1103 memorySpace> &nodalData,
1104 ValueTypeBasisCoeff * quadratureValues,
1105 ValueTypeBasisCoeff * quadratureGradients,
1106 const std::pair<unsigned int, unsigned int> cellRange)
const;
1122 const ValueTypeBasisCoeff * nodalData,
1123 ValueTypeBasisCoeff * quadratureValues,
1124 ValueTypeBasisCoeff * quadratureGradients,
1125 const std::pair<unsigned int, unsigned int> cellRange)
const;
1141 const ValueTypeBasisCoeff *quadratureValues,
1142 const ValueTypeBasisCoeff *quadratureGradients,
1146 & mapQuadIdToProcId,
1147 const std::pair<unsigned int, unsigned int> cellRange)
const;
1162 memorySpace> &nodalData,
1163 ValueTypeBasisCoeff * cellNodalDataPtr,
1164 const std::pair<unsigned int, unsigned int> cellRange)
const;
1177 const ValueTypeBasisCoeff *cellNodalDataPtr,
1180 const std::pair<unsigned int, unsigned int> cellRange)
const;
unsigned int d_quadratureIndex
Definition FEBasisOperations.h:1011
dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > d_flattenedCellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:866
const auto & massVector() const
diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:623
const dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > & quadPoints() const
quad point coordinates for each cell.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseSqrtMassVectorBasisData() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisData.
const auto & cellMassVector() const
Cell level diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:551
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseMassVectorCoeffType
Definition FEBasisOperations.h:952
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:972
std::vector< unsigned int > d_quadratureIDsVector
Definition FEBasisOperations.h:1009
const auto & inverseSqrtMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:569
unsigned int nQuadsPerCell() const
Number of quadrature points per cell for the quadratureID set in reinit.
void computeCellMassMatrix(const unsigned int quadratureID, const unsigned int cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
bool areAllCellsAffine
Definition FEBasisOperations.h:1020
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:956
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseStiffnessVectorBasisData() const
diagonal inverse stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< dftfe::global_size_type, dftfe::utils::MemorySpace::HOST > d_cellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:860
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:993
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_cellInverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:926
unsigned int cellsTypeFlag() const
returns 2 if all cells on current processor are Cartesian, 1 if all cells on current processor are af...
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtMassVectorCoeffType
Definition FEBasisOperations.h:940
const dealii::MatrixFree< 3, ValueTypeBasisData > * d_matrixFreeDataPtr
Definition FEBasisOperations.h:857
const auto & cellInverseMassVector() const
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:517
unsigned int d_nDofsPerCell
Definition FEBasisOperations.h:1017
void initializeMPIPattern()
Constructs the MPIPatternP2P object.
unsigned int d_nVectors
Definition FEBasisOperations.h:1014
void interpolate(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients=NULL) const
Interpolate process level nodal data to cell level quadrature data.
std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > d_cellIndexToCellIteratorMap
Definition FEBasisOperations.h:869
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock
Definition FEBasisOperations.h:90
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseMassVectorBasisType
Definition FEBasisOperations.h:921
void reinit(const unsigned int &vecBlockSize, const unsigned int &cellBlockSize, const unsigned int &quadratureID, const bool isResizeTempStorageForInerpolation=true, const bool isResizeTempStorageForCellMatrices=false)
sets internal variables and optionally resizes internal temp storage for interpolation operations
const dealii::DoFHandler< 3 > & getDofHandler() const
Return the underlying deal.II dofhandler object.
std::map< dealii::CellId, unsigned int > d_cellIdToCellIndexMap
Definition FEBasisOperations.h:870
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:987
dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > & getMultiVectorSinglePrec(const unsigned int vecBlockSize, const unsigned int index=0) const
Gets single precision scratch multivectors.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_massVectorCoeffType
Definition FEBasisOperations.h:944
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellStiffnessMatrixBasisData() const
Cell level stiffness matrix in ValueTypeBasisData.
void computeWeightedCellMassMatrix(const std::pair< unsigned int, unsigned int > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellMassMatrix) const
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellInverseMassVectorBasisData() const
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
void initializeConstraints()
Initializes the constraintMatrixInfo object.
void integrateWithBasisKernel(const ValueTypeBasisCoeff *quadratureValues, const ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > &mapQuadIdToProcId, const std::pair< unsigned int, unsigned int > cellRange) const
Integrate cell level quadrature data times shape functions to process level nodal data.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_JxWData
Definition FEBasisOperations.h:876
~FEBasisOperations()=default
Default Destructor.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseStiffnessVectorBasisType
Definition FEBasisOperations.h:970
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:976
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientData
Definition FEBasisOperations.h:885
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:974
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & inverseJacobians() const
Inverse Jacobian matrices, for cartesian cells returns the diagonal elements of the inverse Jacobian ...
void extractToCellNodalData(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr) const
Get cell level nodal data from process level nodal data.
unsigned int d_quadratureID
Definition FEBasisOperations.h:1010
const auto & inverseJacobiansBasisData() const
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of...
Definition FEBasisOperations.h:411
const auto & cellStiffnessMatrix() const
Cell level stiffness matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:450
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
Definition FEBasisOperations.h:97
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 re...
dftfe::utils::MemoryStorage< dftfe::global_size_type, dftfe::utils::MemorySpace::HOST > & getFlattenedMapsHost()
const auto & shapeFunctionBasisData(bool transpose=false) const
Shape function values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:358
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_stiffnessVectorBasisType
Definition FEBasisOperations.h:982
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock2
Definition FEBasisOperations.h:90
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:954
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > mpiPatternP2P
Definition FEBasisOperations.h:1025
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:995
void computeInverseSqrtMassVector(const bool basisType=true, const bool ceoffType=false)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:985
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:991
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseStiffnessVectorBasisType
Definition FEBasisOperations.h:980
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessVectorBasisType
Definition FEBasisOperations.h:968
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtMassVectorCoeffType
Definition FEBasisOperations.h:960
void accumulateFromCellNodalData(const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData) const
Accumulate cell level nodal data into process level nodal data.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & massVectorBasisData() const
diagonal mass matrix in ValueTypeBasisData
void reinitializeConstraints(std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector)
Reinitializes the constraintMatrixInfo object.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassMatrixCoeffType
Definition FEBasisOperations.h:919
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtMassVectorBasisType
Definition FEBasisOperations.h:958
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & shapeFunctionGradientData(bool transpose=false) const
Shape function gradient values at quadrature points.
std::vector< dealii::CellId > d_cellIndexToCellIdMap
Definition FEBasisOperations.h:867
void computeWeightedCellNjGradNiMatrix(const std::pair< unsigned int, unsigned int > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiMatrix) const
void clearScratchMultiVectors() const
Clears scratch multivectors.
std::vector< UpdateFlags > d_updateFlags
Definition FEBasisOperations.h:1022
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_inverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:950
void interpolateKernel(const ValueTypeBasisCoeff *nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< unsigned int, unsigned int > cellRange) const
Interpolate cell level nodal data to cell level quadrature data.
unsigned int nVectors() const
Number of vectors set in reinit.
const auto & sqrtMassVector() const
sqrt diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:606
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:999
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & JxW() const
determinant of Jacobian times the quadrature weight at each quad point for each cell.
unsigned int nCells() const
Number of locally owned cells on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessMatrixBasisType
Definition FEBasisOperations.h:913
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & cellInverseMassVectorBasisDataSinglePrec() const
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessVectorCoeffType
Definition FEBasisOperations.h:989
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:930
unsigned int cellIndex(const dealii::CellId cellid) const
returns the cell index corresponding to given deal.ii cellID.
void clear()
Clears the FEBasisOperations internal storage.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisDataTranspose
Definition FEBasisOperations.h:910
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > > d_quadPoints
Definition FEBasisOperations.h:864
void distribute(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector, unsigned int constraintIndex=std::numeric_limits< unsigned int >::max()) const
Apply constraints on given multivector.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellSqrtMassVectorBasisData() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseMassVectorBasisType
Definition FEBasisOperations.h:946
const auto & cellSqrtMassVector() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:534
unsigned int d_dofHandlerID
Definition FEBasisOperations.h:1013
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisData
Definition FEBasisOperations.h:904
void computeStiffnessVector(const bool basisType=true, const bool ceoffType=false)
Computes the stiffness matrix \grad Ni . \grad Ni.
unsigned int d_locallyOwnedSize
Definition FEBasisOperations.h:1019
void computeWeightedCellNjGradNiPlusNiGradNjMatrix(const std::pair< unsigned int, unsigned int > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsData
Definition FEBasisOperations.h:87
dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > zeroIndexVec
Definition FEBasisOperations.h:93
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataInternalLayout
Definition FEBasisOperations.h:882
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellValuesBlockCoeff
Definition FEBasisOperations.h:95
void accumulateFromCellNodalDataKernel(const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, const std::pair< unsigned int, unsigned int > cellRange) const
Accumulate cell level nodal data into process level nodal data.
unsigned int nRelaventDofs() const
Number of DoFs on the current processor, locally owned + ghosts.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessMatrixCoeffType
Definition FEBasisOperations.h:915
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionData
Definition FEBasisOperations.h:879
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_stiffnessVectorCoeffType
Definition FEBasisOperations.h:997
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellMassVectorBasisData() const
Cell level diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:932
std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > * d_constraintsVector
Definition FEBasisOperations.h:856
const auto & inverseMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:588
void init(dealii::MatrixFree< 3, ValueTypeBasisData > &matrixFreeData, std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector, const unsigned int &dofHandlerID, const std::vector< unsigned int > &quadratureID, const std::vector< UpdateFlags > updateFlags)
fills required data structures for the given dofHandlerID
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_JxWBasisData
Definition FEBasisOperations.h:898
unsigned int d_cellsBlockSize
Definition FEBasisOperations.h:1016
void computeCellStiffnessMatrix(const unsigned int quadratureID, const unsigned int cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
Computes the cell-level stiffness matrix.
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & inverseMassVectorBasisDataSinglePrec() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellValuesBlock
Definition FEBasisOperations.h:90
dealii::DoFHandler< 3 >::active_cell_iterator getCellIterator(const unsigned int iElem) const
returns the deal.ii cell_iterator corresponing to given cell Index.
void initializeShapeFunctionAndJacobianData()
Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtMassVectorBasisType
Definition FEBasisOperations.h:938
void initializeShapeFunctionAndJacobianBasisData()
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseMassVectorCoeffType
Definition FEBasisOperations.h:928
unsigned int nDofsPerCell() const
Number of DoFs per cell for the dofHandlerID set in init.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:978
const auto & JxWBasisData() const
determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each quad point for each...
Definition FEBasisOperations.h:433
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisData
Definition FEBasisOperations.h:901
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtStiffnessVectorBasisData() const
diagonal sqrt stiffness matrix in ValueTypeBasisData
void init(const FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpaceSrc > &basisOperationsSrc)
fills required data structures from another FEBasisOperations object
std::map< unsigned int, std::vector< dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > > > scratchMultiVectorsSinglePrec
Definition FEBasisOperations.h:1007
bool areAllCellsCartesian
Definition FEBasisOperations.h:1021
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellNodalData
Definition FEBasisOperations.h:87
void createScratchMultiVectors(const unsigned int vecBlockSize, const unsigned int numMultiVecs=1) const
Creates scratch multivectors.
unsigned int d_nCells
Definition FEBasisOperations.h:1015
unsigned int d_localSize
Definition FEBasisOperations.h:1018
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassVectorBasisType
Definition FEBasisOperations.h:934
void initializeFlattenedIndexMaps()
Initializes indexset maps from process level indices to cell level indices for multivectors.
const auto & cellInverseSqrtMassVector() const
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:500
void initializeIndexMaps()
Initializes indexset maps from process level indices to cell level indices for a single vector,...
std::vector< dftUtils::constraintMatrixInfo< memorySpace > > d_constraintInfo
Definition FEBasisOperations.h:853
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseSqrtStiffnessVectorBasisData() const
diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionDataTranspose
Definition FEBasisOperations.h:888
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_inverseJacobianData
Definition FEBasisOperations.h:873
FEBasisOperations(std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr)
Constructor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassVectorCoeffType
Definition FEBasisOperations.h:936
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & getMultiVector(const unsigned int vecBlockSize, const unsigned int index=0) const
Gets scratch multivectors.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseMassVectorBasisData() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellMatrixBlock
Definition FEBasisOperations.h:91
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_inverseJacobianBasisData
Definition FEBasisOperations.h:895
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtMassVectorBasisData() const
sqrt diagonal mass matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassMatrixBasisType
Definition FEBasisOperations.h:917
void createMultiVectorSinglePrec(const unsigned int blocksize, dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > &multiVector) const
Creates a multivector.
void createScratchMultiVectorsSinglePrec(const unsigned int vecBlockSize, const unsigned int numMultiVecs=1) const
Creates single precision scratch multivectors.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & stiffnessVectorBasisData() const
diagonal stiffness matrix in ValueTypeBasisData
dealii::CellId cellID(const unsigned int iElem) const
returns the deal.ii cellID corresponing to given cell Index.
std::map< unsigned int, std::vector< dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > > > scratchMultiVectors
Definition FEBasisOperations.h:965
void createMultiVector(const unsigned int blocksize, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector) const
Creates a multivector.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsDataNonAffine
Definition FEBasisOperations.h:88
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & shapeFunctionData(bool transpose=false) const
Shape function values at quadrature points.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataTranspose
Definition FEBasisOperations.h:891
std::vector< unsigned int > d_nQuadsPerCell
Definition FEBasisOperations.h:1012
const auto & cellMassMatrix() const
Cell level mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:475
const auto & shapeFunctionGradientBasisData(bool transpose=false) const
Shape function gradient values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:383
unsigned int d_nOMPThreads
Definition FEBasisOperations.h:854
const dealii::MatrixFree< 3, ValueTypeBasisData > & matrixFreeData() const
Return the underlying deal.II matrixfree object.
void extractToCellNodalDataKernel(const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr, const std::pair< unsigned int, unsigned int > cellRange) const
Get cell level nodal data from process level nodal data.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellInverseSqrtMassVectorBasisData() const
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellMassMatrixBasisData() const
Cell level mass matrix in ValueTypeBasisData.
void integrateWithBasis(ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > &mapQuadIdToProcId) const
Integrate cell level quadrature data times shape functions to process level nodal data.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_massVectorBasisType
Definition FEBasisOperations.h:942
void interpolateKernel(const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< unsigned int, unsigned int > cellRange) const
Interpolate process level nodal data to cell level quadrature data.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisDataTranspose
Definition FEBasisOperations.h:907
unsigned int nOwnedDofs() const
Number of locally owned DoFs on the current processor.
Definition BLASWrapper.h:35
An class template to encapsulate a MultiVector. A MultiVector is a collection of vectors belonging t...
Definition MultiVector.h:127
Definition MemoryStorage.h:33
Definition FEBasisOperations.h:30
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
Definition FEBasisOperations.h:49
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:58
UpdateFlags
Definition FEBasisOperations.h:32
@ update_values
Definition FEBasisOperations.h:35
@ update_quadpoints
Definition FEBasisOperations.h:41
@ update_transpose
Definition FEBasisOperations.h:39
@ update_inversejacobians
Definition FEBasisOperations.h:43
@ update_jxw
Definition FEBasisOperations.h:45
@ update_default
Definition FEBasisOperations.h:33
@ update_gradients
Definition FEBasisOperations.h:37
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
Definition FEBasisOperations.h:65
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:73
MemorySpace
Definition MemorySpaceType.h:33
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
unsigned long int global_size_type
Definition TypeConfig.h:7
T type
Definition dftfeDataTypes.h:115