18#ifndef dftfeFEBasisOperations_h
19#define dftfeFEBasisOperations_h
83 template <
typename ValueTypeBasisCoeff,
84 typename ValueTypeBasisData,
99 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
136 std::vector<
const dealii::AffineConstraints<ValueTypeBasisData> *>
139 const std::vector<dftfe::uInt> &quadratureID,
140 const std::vector<UpdateFlags> updateFlags);
146 template <dftfe::utils::MemorySpace memorySpaceSrc>
150 memorySpaceSrc> &basisOperationsSrc);
173 const bool isResizeTempStorageForInerpolation =
true,
174 const bool isResizeTempStorageForCellMatrices =
false,
175 const bool isResizeTempStorageForIntegralEvaluations =
false);
207 std::vector<
const dealii::AffineConstraints<ValueTypeBasisData> *>
235 const bool basisType =
false,
236 const bool ceoffType =
true);
241 const bool basisType =
false,
242 const bool ceoffType =
true);
246 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
249 &weightedCellMassMatrix)
const;
253 const std::vector<dftfe::uInt> &cellIndices,
261 &scalarFieldTimesShapeFunctionIntegral)
const;
265 const std::vector<dftfe::uInt> &cellIndices,
273 &scalarFieldTimesGradientShapeFunctionIntegral)
const;
277 const std::vector<dftfe::uInt> &cellIndices,
285 &vectorFieldDyadicGradientShapeFunctionIntegral)
const;
289 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
292 &weightedCellNjGradNiMatrix)
const;
296 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
299 &weightedCellNjGradNiPlusNiGradNjMatrix)
const;
303 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
306 &weightedCellNjGradNiPlusNiGradNjMatrix)
const;
310 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
313 &weightedCellStiffnessMatrix)
const;
317 const bool ceoffType =
false);
325 const bool ceoffType =
false);
332 const bool isResizeTempStorageForCellMatrices,
333 const bool isResizeTempStorageForIntegralEvaluations);
434 if constexpr (std::is_same<ValueTypeBasisCoeff,
435 ValueTypeBasisData>::value)
459 if constexpr (std::is_same<ValueTypeBasisCoeff,
460 ValueTypeBasisData>::value)
484 if constexpr (std::is_same<ValueTypeBasisCoeff,
485 ValueTypeBasisData>::value)
507 if constexpr (std::is_same<ValueTypeBasisCoeff,
508 ValueTypeBasisData>::value)
529 if constexpr (std::is_same<ValueTypeBasisCoeff,
530 ValueTypeBasisData>::value)
546 if constexpr (std::is_same<ValueTypeBasisCoeff,
547 ValueTypeBasisData>::value)
571 if constexpr (std::is_same<ValueTypeBasisCoeff,
572 ValueTypeBasisData>::value)
596 if constexpr (std::is_same<ValueTypeBasisCoeff,
597 ValueTypeBasisData>::value)
613 if constexpr (std::is_same<ValueTypeBasisCoeff,
614 ValueTypeBasisData>::value)
630 if constexpr (std::is_same<ValueTypeBasisCoeff,
631 ValueTypeBasisData>::value)
647 if constexpr (std::is_same<ValueTypeBasisCoeff,
648 ValueTypeBasisData>::value)
665 if constexpr (std::is_same<ValueTypeBasisCoeff,
666 ValueTypeBasisData>::value)
684 if constexpr (std::is_same<ValueTypeBasisCoeff,
685 ValueTypeBasisData>::value)
702 if constexpr (std::is_same<ValueTypeBasisCoeff,
703 ValueTypeBasisData>::value)
719 if constexpr (std::is_same<ValueTypeBasisCoeff,
720 ValueTypeBasisData>::value)
839 dealii::DoFHandler<3>::active_cell_iterator
870 memorySpace> &multiVector)
const;
927 memorySpace> &multiVector,
929 std::numeric_limits<dftfe::uInt>::max())
const;
936 const dealii::MatrixFree<3, ValueTypeBasisData> &
942 const dealii::DoFHandler<3> &
949 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
964 std::vector<dealii::DoFHandler<3>::active_cell_iterator>
1102 std::map<std::pair<dftfe::uInt, dftfe::uInt>,
1127 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
1137 &quadratureValueData,
1139 &quadratureGradValueData,
1141 &quadratureHessianValueData,
1142 const bool isEvaluateGradData =
false,
1143 const bool isEvaluateHessianData =
false,
1144 const bool isEvaluateData =
true)
const;
1152 &quadratureValueData,
1154 &quadratureGradValueData,
1156 &quadratureHessianValueData,
1157 const bool isEvaluateGradData =
false,
1158 const bool isEvaluateHessianData =
false,
1159 const bool isEvaluateData =
true)
const;
1173 memorySpace> &nodalData,
1174 ValueTypeBasisCoeff *quadratureValues,
1175 ValueTypeBasisCoeff *quadratureGradients = NULL)
const;
1189 ValueTypeBasisCoeff *quadratureValues,
1190 ValueTypeBasisCoeff *quadratureGradients,
1194 &mapQuadIdToProcId)
const;
1207 ValueTypeBasisCoeff *cellNodalDataPtr)
const;
1217 const ValueTypeBasisCoeff *cellNodalDataPtr,
1236 memorySpace> &nodalData,
1237 ValueTypeBasisCoeff *quadratureValues,
1238 ValueTypeBasisCoeff *quadratureGradients,
1239 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1255 const ValueTypeBasisCoeff *nodalData,
1256 ValueTypeBasisCoeff *quadratureValues,
1257 ValueTypeBasisCoeff *quadratureGradients,
1258 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1282 const ValueTypeBasisCoeff *quadratureValues,
1283 const ValueTypeBasisCoeff *quadratureGradients,
1288 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1303 memorySpace> &nodalData,
1304 ValueTypeBasisCoeff *cellNodalDataPtr,
1305 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
1318 const ValueTypeBasisCoeff *cellNodalDataPtr,
1321 const std::pair<dftfe::uInt, dftfe::uInt> cellRange)
const;
dftfe::uInt d_quadratureIndex
Definition FEBasisOperations.h:1114
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:1112
void shapeFunctionsCenteredAtQuad1EvaluatedAtQuad2(const dftfe::uInt quadId1, const dftfe::uInt quadId2)
const auto & massVector() const
diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:717
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:645
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseMassVectorCoeffType
Definition FEBasisOperations.h:1054
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:1074
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisData
Definition FEBasisOperations.h:1000
const auto & inverseSqrtMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:663
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & collocationShapeFunctionGradientData() const
Collocation shape function gradient values at quadrature points.
bool areAllCellsAffine
Definition FEBasisOperations.h:1123
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:1058
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseStiffnessVectorBasisData() const
diagonal inverse stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1095
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_cellInverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:1028
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtMassVectorCoeffType
Definition FEBasisOperations.h:1042
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_inverseJacobianData
Definition FEBasisOperations.h:969
const dealii::MatrixFree< 3, ValueTypeBasisData > * d_matrixFreeDataPtr
Definition FEBasisOperations.h:951
const auto & cellInverseMassVector() const
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:611
void initializeMPIPattern()
Constructs the MPIPatternP2P object.
dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > d_cellCentroids
Definition FEBasisOperations.h:960
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:965
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock
Definition FEBasisOperations.h:93
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseMassVectorBasisType
Definition FEBasisOperations.h:1023
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:1122
std::map< dealii::CellId, dftfe::uInt > d_cellIdToCellIndexMap
Definition FEBasisOperations.h:966
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > > d_quadPoints
Definition FEBasisOperations.h:957
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:1089
dftfe::uInt nOwnedDofs() const
Number of locally owned DoFs on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_massVectorCoeffType
Definition FEBasisOperations.h:1046
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataInternalLayout
Definition FEBasisOperations.h:978
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:990
~FEBasisOperations()=default
Default Destructor.
dftfe::uInt d_dofHandlerID
Definition FEBasisOperations.h:1116
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseStiffnessVectorBasisType
Definition FEBasisOperations.h:1072
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionData
Definition FEBasisOperations.h:975
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:1078
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:1076
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:987
const auto & inverseJacobiansBasisData() const
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of...
Definition FEBasisOperations.h:505
const auto & cellStiffnessMatrix() const
Cell level stiffness matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:544
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
Definition FEBasisOperations.h:100
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:1118
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > zeroIndexVec
Definition FEBasisOperations.h:96
const auto & shapeFunctionBasisData(bool transpose=false) const
Shape function values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:432
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_stiffnessVectorBasisType
Definition FEBasisOperations.h:1084
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock2
Definition FEBasisOperations.h:93
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:1056
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > mpiPatternP2P
Definition FEBasisOperations.h:1128
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1097
void computeInverseSqrtMassVector(const bool basisType=true, const bool ceoffType=false)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:1087
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:1093
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseStiffnessVectorBasisType
Definition FEBasisOperations.h:1082
dftfe::uInt d_cellsBlockSize
Definition FEBasisOperations.h:1119
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessVectorBasisType
Definition FEBasisOperations.h:1070
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtMassVectorCoeffType
Definition FEBasisOperations.h:1062
dftfe::uInt d_nVectors
Definition FEBasisOperations.h:1117
dftfe::uInt d_nOMPThreads
Definition FEBasisOperations.h:948
std::vector< dftfe::uInt > d_nQuadsPerCell
Definition FEBasisOperations.h:1115
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.
void computeScalarFieldTimesGradientShapeFunctionIntegral(const std::vector< dftfe::uInt > &cellIndices, const dftfe::uInt &noKpoints, const dftfe::uInt &noOfVectors, const dftfe::uInt &totalElements, const dftfe::uInt &iElemStart, const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > &scalarField, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > &scalarFieldTimesGradientShapeFunctionIntegral) const
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassMatrixCoeffType
Definition FEBasisOperations.h:1021
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:1060
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:963
void clearScratchMultiVectors() const
Clears scratch multivectors.
std::vector< UpdateFlags > d_updateFlags
Definition FEBasisOperations.h:1125
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_inverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:1052
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_collocationShapeFunctionGradientBasisData
Definition FEBasisOperations.h:1003
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 computeScalarFieldTimesShapeFunctionIntegral(const std::vector< dftfe::uInt > &cellIndices, const dftfe::uInt &noKpoints, const dftfe::uInt &noOfVectors, const dftfe::uInt &totalElements, const dftfe::uInt &iElemStart, const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > &scalarField, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > &scalarFieldTimesShapeFunctionIntegral) const
const auto & sqrtMassVector() const
sqrt diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:700
dftfe::uInt d_quadratureID
Definition FEBasisOperations.h:1113
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_JxWData
Definition FEBasisOperations.h:972
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:1101
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:1110
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_collocationShapeFunctionGradientData
Definition FEBasisOperations.h:981
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessMatrixBasisType
Definition FEBasisOperations.h:1015
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > d_flattenedCellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:962
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:1091
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:1032
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:984
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:1048
const auto & cellSqrtMassVector() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:628
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:1009
void computeStiffnessVector(const bool basisType=true, const bool ceoffType=false)
Computes the stiffness matrix \grad Ni . \grad Ni.
void computeVectorFieldDyadicGradientShapeFunctionIntegral(const std::vector< dftfe::uInt > &cellIndices, const dftfe::uInt &noKpoints, const dftfe::uInt &noOfVectors, const dftfe::uInt &totalElements, const dftfe::uInt &iElemStart, const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > &vectorField, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > &vectorFieldDyadicGradientShapeFunctionIntegral) const
void interpolateNoConstraints(const distributedCPUVec< double > &nodalField, const dftfe::uInt dofHandlerId, const dftfe::uInt quadratureId, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureGradValueData, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureHessianValueData, const bool isEvaluateGradData=false, const bool isEvaluateHessianData=false, const bool isEvaluateData=true) const
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:90
std::map< dftfe::uInt, std::vector< dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > > > scratchMultiVectors
Definition FEBasisOperations.h:1067
dftfe::uInt nCells() const
Number of locally owned cells on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellValuesBlockCoeff
Definition FEBasisOperations.h:98
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessMatrixCoeffType
Definition FEBasisOperations.h:1017
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_stiffnessVectorCoeffType
Definition FEBasisOperations.h:1099
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:1034
std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > * d_constraintsVector
Definition FEBasisOperations.h:950
const auto & inverseMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:682
dftfe::uInt d_localSize
Definition FEBasisOperations.h:1121
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:93
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:1040
void initializeShapeFunctionAndJacobianBasisData()
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseMassVectorCoeffType
Definition FEBasisOperations.h:1030
void reinit(const dftfe::uInt &vecBlockSize, const dftfe::uInt &cellBlockSize, const dftfe::uInt &quadratureID, const bool isResizeTempStorageForInerpolation=true, const bool isResizeTempStorageForCellMatrices=false, const bool isResizeTempStorageForIntegralEvaluations=false)
sets internal variables and optionally resizes internal temp storage for interpolation operations
std::map< std::pair< dftfe::uInt, dftfe::uInt >, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFnValQuad1ToQuad2
Definition FEBasisOperations.h:1104
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:1080
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:527
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:1124
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellNodalData
Definition FEBasisOperations.h:90
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassVectorBasisType
Definition FEBasisOperations.h:1036
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:594
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:953
dftfe::uInt d_nDofsPerCell
Definition FEBasisOperations.h:1120
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellGradientsBlockCoeff
Definition FEBasisOperations.h:98
std::vector< dftUtils::constraintMatrixInfo< memorySpace > > d_constraintInfo
Definition FEBasisOperations.h:947
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:1006
FEBasisOperations(std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr)
Constructor.
void interpolate(distributedCPUVec< double > &nodalField, const dftfe::uInt dofHandlerId, const dftfe::uInt quadratureId, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureGradValueData, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureHessianValueData, const bool isEvaluateGradData=false, const bool isEvaluateHessianData=false, const bool isEvaluateData=true) const
const auto & collocationShapeFunctionGradientBasisData() const
Collocation shape function gradient values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:482
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:1038
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_JxWBasisData
Definition FEBasisOperations.h:997
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseMassVectorBasisData() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellMatrixBlock
Definition FEBasisOperations.h:94
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtMassVectorBasisData() const
sqrt diagonal mass matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassMatrixBasisType
Definition FEBasisOperations.h:1019
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
const dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > & cellCentroids() const
quad point coordinates for each cell.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsDataNonAffine
Definition FEBasisOperations.h:91
void interpolateQ1ToQ2(const double *Q1Field, const dftfe::uInt quadRule1, const dftfe::uInt quadRule2, double *Q2Field, const dftfe::uInt numComponents) const
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisDataTranspose
Definition FEBasisOperations.h:1012
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:569
const auto & shapeFunctionGradientBasisData(bool transpose=false) const
Shape function gradient values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:457
const dealii::MatrixFree< 3, ValueTypeBasisData > & matrixFreeData() const
Return the underlying deal.II matrixfree object.
void resizeTempStorage(const bool isResizeTempStorageForInerpolation, const bool isResizeTempStorageForCellMatrices, const bool isResizeTempStorageForIntegralEvaluations)
Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in re...
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:1044
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_inverseJacobianBasisData
Definition FEBasisOperations.h:994
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:51
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:60
UpdateFlags
Definition FEBasisOperations.h:32
@ update_values
Definition FEBasisOperations.h:35
@ update_quadpoints
Definition FEBasisOperations.h:41
@ update_collocation_gradients
Definition FEBasisOperations.h:47
@ 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:68
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:76
MemorySpace
Definition MemorySpaceType.h:33
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
dealii::LinearAlgebra::distributed::Vector< elem_type, dealii::MemorySpace::Host > distributedCPUVec
Definition headers.h:92
T type
Definition dftfeDataTypes.h:112