18#ifndef dftfeFEBasisOperations_h
19#define dftfeFEBasisOperations_h
81 template <
typename ValueTypeBasisCoeff,
82 typename ValueTypeBasisData,
97 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
134 std::vector<
const dealii::AffineConstraints<ValueTypeBasisData> *>
137 const std::vector<dftfe::uInt> &quadratureID,
138 const std::vector<UpdateFlags> updateFlags);
144 template <dftfe::utils::MemorySpace memorySpaceSrc>
148 memorySpaceSrc> &basisOperationsSrc);
165 const bool isResizeTempStorageForInerpolation =
true,
166 const bool isResizeTempStorageForCellMatrices =
false);
198 std::vector<
const dealii::AffineConstraints<ValueTypeBasisData> *>
226 const bool basisType =
false,
227 const bool ceoffType =
true);
232 const bool basisType =
false,
233 const bool ceoffType =
true);
237 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
240 &weightedCellMassMatrix)
const;
244 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
247 &weightedCellNjGradNiMatrix)
const;
251 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
254 &weightedCellNjGradNiPlusNiGradNjMatrix)
const;
258 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
261 &weightedCellNjGradNiPlusNiGradNjMatrix)
const;
265 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
268 &weightedCellStiffnessMatrix)
const;
272 const bool ceoffType =
false);
280 const bool ceoffType =
false);
287 const bool isResizeTempStorageForCellMatrices);
375 if constexpr (std::is_same<ValueTypeBasisCoeff,
376 ValueTypeBasisData>::value)
400 if constexpr (std::is_same<ValueTypeBasisCoeff,
401 ValueTypeBasisData>::value)
428 if constexpr (std::is_same<ValueTypeBasisCoeff,
429 ValueTypeBasisData>::value)
450 if constexpr (std::is_same<ValueTypeBasisCoeff,
451 ValueTypeBasisData>::value)
467 if constexpr (std::is_same<ValueTypeBasisCoeff,
468 ValueTypeBasisData>::value)
492 if constexpr (std::is_same<ValueTypeBasisCoeff,
493 ValueTypeBasisData>::value)
517 if constexpr (std::is_same<ValueTypeBasisCoeff,
518 ValueTypeBasisData>::value)
534 if constexpr (std::is_same<ValueTypeBasisCoeff,
535 ValueTypeBasisData>::value)
551 if constexpr (std::is_same<ValueTypeBasisCoeff,
552 ValueTypeBasisData>::value)
568 if constexpr (std::is_same<ValueTypeBasisCoeff,
569 ValueTypeBasisData>::value)
586 if constexpr (std::is_same<ValueTypeBasisCoeff,
587 ValueTypeBasisData>::value)
605 if constexpr (std::is_same<ValueTypeBasisCoeff,
606 ValueTypeBasisData>::value)
623 if constexpr (std::is_same<ValueTypeBasisCoeff,
624 ValueTypeBasisData>::value)
640 if constexpr (std::is_same<ValueTypeBasisCoeff,
641 ValueTypeBasisData>::value)
760 dealii::DoFHandler<3>::active_cell_iterator
791 memorySpace> &multiVector)
const;
848 memorySpace> &multiVector,
850 std::numeric_limits<dftfe::uInt>::max())
const;
857 const dealii::MatrixFree<3, ValueTypeBasisData> &
863 const dealii::DoFHandler<3> &
870 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
882 std::vector<dealii::DoFHandler<3>::active_cell_iterator>
1038 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
1054 memorySpace> &nodalData,
1055 ValueTypeBasisCoeff *quadratureValues,
1056 ValueTypeBasisCoeff *quadratureGradients = NULL)
const;
1070 ValueTypeBasisCoeff *quadratureValues,
1071 ValueTypeBasisCoeff *quadratureGradients,
1075 &mapQuadIdToProcId)
const;
1088 ValueTypeBasisCoeff *cellNodalDataPtr)
const;
1098 const ValueTypeBasisCoeff *cellNodalDataPtr,
1117 memorySpace> &nodalData,
1118 ValueTypeBasisCoeff *quadratureValues,
1119 ValueTypeBasisCoeff *quadratureGradients,
1120 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1136 const ValueTypeBasisCoeff *nodalData,
1137 ValueTypeBasisCoeff *quadratureValues,
1138 ValueTypeBasisCoeff *quadratureGradients,
1139 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1155 const ValueTypeBasisCoeff *quadratureValues,
1156 const ValueTypeBasisCoeff *quadratureGradients,
1161 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1176 memorySpace> &nodalData,
1177 ValueTypeBasisCoeff *cellNodalDataPtr,
1178 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1191 const ValueTypeBasisCoeff *cellNodalDataPtr,
1194 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
dftfe::uInt d_quadratureIndex
Definition FEBasisOperations.h:1025
dftfe::uInt nVectors() const
Number of vectors set in reinit.
void computeWeightedCellStiffnessMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellStiffnessMatrix) const
std::vector< dftfe::uInt > d_quadratureIDsVector
Definition FEBasisOperations.h:1023
const auto & massVector() const
diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:638
dealii::CellId cellID(const dftfe::uInt iElem) const
returns the deal.ii cellID corresponing to given cell Index.
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:566
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseMassVectorCoeffType
Definition FEBasisOperations.h:966
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:986
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisData
Definition FEBasisOperations.h:915
const auto & inverseSqrtMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:584
bool areAllCellsAffine
Definition FEBasisOperations.h:1034
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:970
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseStiffnessVectorBasisData() const
diagonal inverse stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1007
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_cellInverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:940
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtMassVectorCoeffType
Definition FEBasisOperations.h:954
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_inverseJacobianData
Definition FEBasisOperations.h:887
const dealii::MatrixFree< 3, ValueTypeBasisData > * d_matrixFreeDataPtr
Definition FEBasisOperations.h:872
const auto & cellInverseMassVector() const
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:532
void initializeMPIPattern()
Constructs the MPIPatternP2P 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.
std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > d_cellIndexToCellIteratorMap
Definition FEBasisOperations.h:883
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock
Definition FEBasisOperations.h:91
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseMassVectorBasisType
Definition FEBasisOperations.h:935
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.
dftfe::uInt d_locallyOwnedSize
Definition FEBasisOperations.h:1033
std::map< dealii::CellId, dftfe::uInt > d_cellIdToCellIndexMap
Definition FEBasisOperations.h:884
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > > d_quadPoints
Definition FEBasisOperations.h:878
const dealii::DoFHandler< 3 > & getDofHandler() const
Return the underlying deal.II dofhandler object.
void computeCellMassMatrix(const dftfe::uInt quadratureID, const dftfe::uInt cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1001
dftfe::uInt nOwnedDofs() const
Number of locally owned DoFs on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_massVectorCoeffType
Definition FEBasisOperations.h:958
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataInternalLayout
Definition FEBasisOperations.h:896
dealii::DoFHandler< 3 >::active_cell_iterator getCellIterator(const dftfe::uInt iElem) const
returns the deal.ii cell_iterator corresponing to given cell Index.
void computeWeightedCellNjGradNiMinusNiGradNjMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellStiffnessMatrixBasisData() const
Cell level stiffness matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellInverseMassVectorBasisData() const
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
void initializeConstraints()
Initializes the constraintMatrixInfo object.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataTranspose
Definition FEBasisOperations.h:905
~FEBasisOperations()=default
Default Destructor.
dftfe::uInt d_dofHandlerID
Definition FEBasisOperations.h:1027
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseStiffnessVectorBasisType
Definition FEBasisOperations.h:984
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionData
Definition FEBasisOperations.h:893
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:990
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:988
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.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionDataTranspose
Definition FEBasisOperations.h:902
const auto & inverseJacobiansBasisData() const
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of...
Definition FEBasisOperations.h:426
const auto & cellStiffnessMatrix() const
Cell level stiffness matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:465
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
Definition FEBasisOperations.h:98
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...
void computeCellStiffnessMatrix(const dftfe::uInt quadratureID, const dftfe::uInt cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
Computes the cell-level stiffness matrix.
dftfe::uInt d_nCells
Definition FEBasisOperations.h:1029
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > zeroIndexVec
Definition FEBasisOperations.h:94
const auto & shapeFunctionBasisData(bool transpose=false) const
Shape function values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:373
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< ValueTypeBasisData, memorySpace > d_stiffnessVectorBasisType
Definition FEBasisOperations.h:996
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock2
Definition FEBasisOperations.h:91
void createMultiVector(const dftfe::uInt blocksize, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector) const
Creates a multivector.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:968
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > mpiPatternP2P
Definition FEBasisOperations.h:1039
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1009
void computeInverseSqrtMassVector(const bool basisType=true, const bool ceoffType=false)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:999
void createMultiVectorSinglePrec(const dftfe::uInt blocksize, dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > &multiVector) const
Creates a multivector.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1005
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseStiffnessVectorBasisType
Definition FEBasisOperations.h:994
dftfe::uInt d_cellsBlockSize
Definition FEBasisOperations.h:1030
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessVectorBasisType
Definition FEBasisOperations.h:982
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtMassVectorCoeffType
Definition FEBasisOperations.h:974
dftfe::uInt d_nVectors
Definition FEBasisOperations.h:1028
dftfe::uInt d_nOMPThreads
Definition FEBasisOperations.h:869
std::vector< dftfe::uInt > d_nQuadsPerCell
Definition FEBasisOperations.h:1026
dftfe::uInt nDofsPerCell() const
Number of DoFs per cell for the dofHandlerID set in init.
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:933
void computeWeightedCellNjGradNiPlusNiGradNjMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtMassVectorBasisType
Definition FEBasisOperations.h:972
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:881
void clearScratchMultiVectors() const
Clears scratch multivectors.
std::vector< UpdateFlags > d_updateFlags
Definition FEBasisOperations.h:1036
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_inverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:964
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.
const auto & sqrtMassVector() const
sqrt diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:621
dftfe::uInt d_quadratureID
Definition FEBasisOperations.h:1024
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_JxWData
Definition FEBasisOperations.h:890
void createScratchMultiVectors(const dftfe::uInt vecBlockSize, const dftfe::uInt numMultiVecs=1) const
Creates scratch multivectors.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:1013
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & JxW() const
determinant of Jacobian times the quadrature weight at each quad point for each cell.
std::map< dftfe::uInt, std::vector< dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > > > scratchMultiVectorsSinglePrec
Definition FEBasisOperations.h:1021
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessMatrixBasisType
Definition FEBasisOperations.h:927
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > d_flattenedCellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:880
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:1003
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:944
void clear()
Clears the FEBasisOperations internal storage.
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 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.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientData
Definition FEBasisOperations.h:899
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellSqrtMassVectorBasisData() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisData.
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.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseMassVectorBasisType
Definition FEBasisOperations.h:960
const auto & cellSqrtMassVector() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:549
void createScratchMultiVectorsSinglePrec(const dftfe::uInt vecBlockSize, const dftfe::uInt numMultiVecs=1) const
Creates single precision scratch multivectors.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisDataTranspose
Definition FEBasisOperations.h:921
void computeStiffnessVector(const bool basisType=true, const bool ceoffType=false)
Computes the stiffness matrix \grad Ni . \grad Ni.
dftfe::uInt cellIndex(const dealii::CellId cellid) const
returns the cell index corresponding to given deal.ii cellID.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsData
Definition FEBasisOperations.h:88
std::map< dftfe::uInt, std::vector< dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > > > scratchMultiVectors
Definition FEBasisOperations.h:979
dftfe::uInt nCells() const
Number of locally owned cells on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellValuesBlockCoeff
Definition FEBasisOperations.h:96
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessMatrixCoeffType
Definition FEBasisOperations.h:929
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_stiffnessVectorCoeffType
Definition FEBasisOperations.h:1011
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:946
std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > * d_constraintsVector
Definition FEBasisOperations.h:871
const auto & inverseMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:603
dftfe::uInt d_localSize
Definition FEBasisOperations.h:1032
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:91
void computeWeightedCellMassMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellMassMatrix) const
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & getMultiVector(const dftfe::uInt vecBlockSize, const dftfe::uInt index=0) const
Gets scratch multivectors.
void initializeShapeFunctionAndJacobianData()
Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtMassVectorBasisType
Definition FEBasisOperations.h:952
void initializeShapeFunctionAndJacobianBasisData()
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseMassVectorCoeffType
Definition FEBasisOperations.h:942
dftfe::uInt nQuadsPerCell() const
Number of quadrature points per cell for the quadratureID set in reinit.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:992
dftfe::uInt cellsTypeFlag() const
returns 2 if all cells on current processor are Cartesian, 1 if all cells on current processor are af...
const auto & JxWBasisData() const
determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each quad point for each...
Definition FEBasisOperations.h:448
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.
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
void computeWeightedCellNjGradNiMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiMatrix) const
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
bool areAllCellsCartesian
Definition FEBasisOperations.h:1035
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellNodalData
Definition FEBasisOperations.h:88
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassVectorBasisType
Definition FEBasisOperations.h:948
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:515
void initializeIndexMaps()
Initializes indexset maps from process level indices to cell level indices for a single vector,...
dftfe::utils::MemoryStorage< dftfe::uInt, dftfe::utils::MemorySpace::HOST > d_cellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:874
dftfe::uInt d_nDofsPerCell
Definition FEBasisOperations.h:1031
std::vector< dftUtils::constraintMatrixInfo< memorySpace > > d_constraintInfo
Definition FEBasisOperations.h:868
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseSqrtStiffnessVectorBasisData() const
diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisData
Definition FEBasisOperations.h:918
FEBasisOperations(std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr)
Constructor.
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.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassVectorCoeffType
Definition FEBasisOperations.h:950
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_JxWBasisData
Definition FEBasisOperations.h:912
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseMassVectorBasisData() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellMatrixBlock
Definition FEBasisOperations.h:92
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtMassVectorBasisData() const
sqrt diagonal mass matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassMatrixBasisType
Definition FEBasisOperations.h:931
void distribute(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector, dftfe::uInt constraintIndex=std::numeric_limits< dftfe::uInt >::max()) const
Apply constraints on given multivector.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & stiffnessVectorBasisData() const
diagonal stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsDataNonAffine
Definition FEBasisOperations.h:89
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisDataTranspose
Definition FEBasisOperations.h:924
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & shapeFunctionData(bool transpose=false) const
Shape function values at quadrature points.
dftfe::uInt nRelaventDofs() const
Number of DoFs on the current processor, locally owned + ghosts.
const auto & cellMassMatrix() const
Cell level mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:490
const auto & shapeFunctionGradientBasisData(bool transpose=false) const
Shape function gradient values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:398
const dealii::MatrixFree< 3, ValueTypeBasisData > & matrixFreeData() const
Return the underlying deal.II matrixfree object.
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.
dftfe::utils::MemoryStorage< dftfe::uInt, dftfe::utils::MemorySpace::HOST > & getFlattenedMapsHost()
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_massVectorBasisType
Definition FEBasisOperations.h:956
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_inverseJacobianBasisData
Definition FEBasisOperations.h:909
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:66
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:74
MemorySpace
Definition MemorySpaceType.h:33
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
T type
Definition dftfeDataTypes.h:113