DFT-FE 1.3.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
FEBasisOperations.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17
18#ifndef dftfeFEBasisOperations_h
19#define dftfeFEBasisOperations_h
20
21#include <MultiVector.h>
22#include <headers.h>
24#include <DeviceTypeConfig.h>
25#include <BLASWrapper.h>
26
27namespace dftfe
28{
29 namespace basis
30 {
32 {
34
35 update_values = 0x0001,
36
38
40
42
44
45 update_jxw = 0x0020,
46
48 };
49
50 inline UpdateFlags
51 operator|(const UpdateFlags f1, const UpdateFlags f2)
52 {
53 return static_cast<UpdateFlags>(static_cast<dftfe::uInt>(f1) |
54 static_cast<dftfe::uInt>(f2));
55 }
56
57
58
59 inline UpdateFlags &
61 {
62 f1 = f1 | f2;
63 return f1;
64 }
65
66
67 inline UpdateFlags
68 operator&(const UpdateFlags f1, const UpdateFlags f2)
69 {
70 return static_cast<UpdateFlags>(static_cast<dftfe::uInt>(f1) &
71 static_cast<dftfe::uInt>(f2));
72 }
73
74
75 inline UpdateFlags &
77 {
78 f1 = f1 & f2;
79 return f1;
80 }
81
82
83 template <typename ValueTypeBasisCoeff,
84 typename ValueTypeBasisData,
85 dftfe::utils::MemorySpace memorySpace>
87 {
88 protected:
99 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
101
102 public:
103 /**
104 * @brief Constructor
105 */
108 BLASWrapperPtr);
109
110
111 /**
112 * @brief Default Destructor
113 */
115
116 /**
117 * @brief Clears the FEBasisOperations internal storage.
118 */
119 void
121
122 /**
123 * @brief fills required data structures for the given dofHandlerID
124 * @param[in] matrixFreeData MatrixFree object.
125 * @param[in] constraintsVector std::vector of AffineConstraints, should
126 * be the same vector which was passed for the construction of the given
127 * MatrixFree object.
128 * @param[in] dofHandlerID dofHandler index to be used for getting data
129 * from the MatrixFree object.
130 * @param[in] quadratureID std::vector of quadratureIDs to be used, should
131 * be the same IDs which were used during the construction of the given
132 * MatrixFree object.
133 */
134 void
135 init(dealii::MatrixFree<3, ValueTypeBasisData> &matrixFreeData,
136 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
137 &constraintsVector,
138 const dftfe::uInt &dofHandlerID,
139 const std::vector<dftfe::uInt> &quadratureID,
140 const std::vector<UpdateFlags> updateFlags);
141
142 /**
143 * @brief fills required data structures from another FEBasisOperations object
144 * @param[in] basisOperationsSrc Source FEBasisOperations object.
145 */
146 template <dftfe::utils::MemorySpace memorySpaceSrc>
147 void
148 init(const FEBasisOperations<ValueTypeBasisCoeff,
149 ValueTypeBasisData,
150 memorySpaceSrc> &basisOperationsSrc);
151
152 /*The shape functions centered at the quadrature points derived from
153 * quadId1 evaluated at quadId2 */
154 void
156 const dftfe::uInt quadId2);
157
158 /**
159 * @brief sets internal variables and optionally resizes internal temp storage for interpolation operations
160 * @param[in] vecBlockSize block size to used for operations on vectors,
161 * this has to be set to the exact value before any such operations are
162 * called.
163 * @param[in] cellBlockSize block size to used for cells, this has to be
164 * set to a value greater than or equal to the required value before any
165 * such operations are called
166 * @param[in] quadratureID Quadrature index to be used.
167 * @param[in] isResizeTempStorage whether to resize internal tempstorage.
168 */
169 void
170 reinit(const dftfe::uInt &vecBlockSize,
171 const dftfe::uInt &cellBlockSize,
172 const dftfe::uInt &quadratureID,
173 const bool isResizeTempStorageForInerpolation = true,
174 const bool isResizeTempStorageForCellMatrices = false,
175 const bool isResizeTempStorageForIntegralEvaluations = false);
176
180
181 // private:
182
183
184 /**
185 * @brief Initializes indexset maps from process level indices to cell level indices for a single vector, also initializes cell index to cellid map.
186 */
187 void
189
190 /**
191 * @brief Initializes indexset maps from process level indices to cell level indices for multivectors.
192 */
193 void
195
196 /**
197 * @brief Initializes the constraintMatrixInfo object.
198 */
199 void
201
202 /**
203 * @brief Reinitializes the constraintMatrixInfo object.
204 */
205 void
207 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
208 &constraintsVector);
209
210 /**
211 * @brief Constructs the MPIPatternP2P object.
212 */
213 void
215
216 /**
217 * @brief Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
218 */
219 void
221
222 /**
223 * @brief Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
224 */
225 void
227
228
229 /**
230 * @brief Computes the cell-level stiffness matrix.
231 */
232 void
234 const dftfe::uInt cellsBlockSize,
235 const bool basisType = false,
236 const bool ceoffType = true);
237
238 void
240 const dftfe::uInt cellsBlockSize,
241 const bool basisType = false,
242 const bool ceoffType = true);
243
244 void
246 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
249 &weightedCellMassMatrix) const;
250
251 void
253 const std::vector<dftfe::uInt> &cellIndices,
254 const dftfe::uInt &noKpoints,
255 const dftfe::uInt &noOfVectors,
256 const dftfe::uInt &totalElements,
257 const dftfe::uInt &iElemStart,
259 &scalarField,
261 &scalarFieldTimesShapeFunctionIntegral) const;
262
263 void
265 const std::vector<dftfe::uInt> &cellIndices,
266 const dftfe::uInt &noKpoints,
267 const dftfe::uInt &noOfVectors,
268 const dftfe::uInt &totalElements,
269 const dftfe::uInt &iElemStart,
271 &scalarField,
273 &scalarFieldTimesGradientShapeFunctionIntegral) const;
274
275 void
277 const std::vector<dftfe::uInt> &cellIndices,
278 const dftfe::uInt &noKpoints,
279 const dftfe::uInt &noOfVectors,
280 const dftfe::uInt &totalElements,
281 const dftfe::uInt &iElemStart,
283 &vectorField,
285 &vectorFieldDyadicGradientShapeFunctionIntegral) const;
286
287 void
289 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
292 &weightedCellNjGradNiMatrix) const;
293
294 void
296 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
299 &weightedCellNjGradNiPlusNiGradNjMatrix) const;
300
301 void
303 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
306 &weightedCellNjGradNiPlusNiGradNjMatrix) const;
307
308 void
310 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
313 &weightedCellStiffnessMatrix) const;
314
315 void
316 computeInverseSqrtMassVector(const bool basisType = true,
317 const bool ceoffType = false);
318
319 /**
320 * @brief Computes the stiffness matrix
321 * \grad Ni . \grad Ni.
322 */
323 void
324 computeStiffnessVector(const bool basisType = true,
325 const bool ceoffType = false);
326
327 /**
328 * @brief Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in reinit.
329 */
330 void
331 resizeTempStorage(const bool isResizeTempStorageForInerpolation,
332 const bool isResizeTempStorageForCellMatrices,
333 const bool isResizeTempStorageForIntegralEvaluations);
334
335 /**
336 * @brief Number of quadrature points per cell for the quadratureID set in reinit.
337 */
340
341 /**
342 * @brief Number of vectors set in reinit.
343 */
345 nVectors() const;
346 /**
347 * @brief Number of DoFs per cell for the dofHandlerID set in init.
348 */
351
352 /**
353 * @brief Number of locally owned cells on the current processor.
354 */
356 nCells() const;
357
358 /**
359 * @brief Number of DoFs on the current processor, locally owned + ghosts.
360 */
363
364 /**
365 * @brief Number of locally owned DoFs on the current processor.
366 */
368 nOwnedDofs() const;
369
370 /**
371 * @brief Shape function values at quadrature points.
372 * @param[in] transpose if false the the data is indexed as [iQuad *
373 * d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
374 * d_nQuadsPerCell + iQuad].
375 */
377 shapeFunctionData(bool transpose = false) const;
378
379 /**
380 * @brief Shape function gradient values at quadrature points.
381 * @param[in] transpose if false the the data is indexed as [iDim *
382 * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
383 * if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
384 * iNode * d_nQuadsPerCell + iQuad].
385 */
387 shapeFunctionGradientData(bool transpose = false) const;
388
389 /**
390 * @brief Collocation shape function gradient values at quadrature points.
391 */
394
395 /**
396 * @brief Inverse Jacobian matrices, for cartesian cells returns the
397 * diagonal elements of the inverse Jacobian matrices for each cell, for
398 * affine cells returns the 3x3 inverse Jacobians for each cell otherwise
399 * returns the 3x3 inverse Jacobians at each quad point for each cell.
400 */
403
404 /**
405 * @brief determinant of Jacobian times the quadrature weight at each
406 * quad point for each cell.
407 */
409 JxW() const;
410
411 /**
412 * @brief quad point coordinates for each cell.
413 */
414 const dftfe::utils::MemoryStorage<ValueTypeBasisData,
416 quadPoints() const;
417
418 /**
419 * @brief quad point coordinates for each cell.
420 */
421 const dftfe::utils::MemoryStorage<ValueTypeBasisData,
424
425 /**
426 * @brief Shape function values at quadrature points in ValueTypeBasisData.
427 * @param[in] transpose if false the the data is indexed as [iQuad *
428 * d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
429 * d_nQuadsPerCell + iQuad].
430 */
431 const auto &
432 shapeFunctionBasisData(bool transpose = false) const
433 {
434 if constexpr (std::is_same<ValueTypeBasisCoeff,
435 ValueTypeBasisData>::value)
436 {
437 return transpose ?
440 }
441 else
442 {
443 return transpose ?
445 ->second :
447 }
448 }
449 /**
450 * @brief Shape function gradient values at quadrature points in ValueTypeBasisData.
451 * @param[in] transpose if false the the data is indexed as [iDim *
452 * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
453 * if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
454 * iNode * d_nQuadsPerCell + iQuad].
455 */
456 const auto &
457 shapeFunctionGradientBasisData(bool transpose = false) const
458 {
459 if constexpr (std::is_same<ValueTypeBasisCoeff,
460 ValueTypeBasisData>::value)
461 {
462 return transpose ?
464 ->second :
466 }
467 else
468 {
469 return transpose ?
471 .find(d_quadratureID)
472 ->second :
474 ->second;
475 }
476 }
477
478 /**
479 * @brief Collocation shape function gradient values at quadrature points in ValueTypeBasisData.
480 */
481 const auto &
483 {
484 if constexpr (std::is_same<ValueTypeBasisCoeff,
485 ValueTypeBasisData>::value)
486 {
488 ->second;
489 }
490 else
491 {
493 .find(d_quadratureID)
494 ->second;
495 }
496 }
497
498 /**
499 * @brief Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the
500 * diagonal elements of the inverse Jacobian matrices for each cell, for
501 * affine cells returns the 3x3 inverse Jacobians for each cell otherwise
502 * returns the 3x3 inverse Jacobians at each quad point for each cell.
503 */
504 const auto &
506 {
507 if constexpr (std::is_same<ValueTypeBasisCoeff,
508 ValueTypeBasisData>::value)
509 {
512 ->second;
513 }
514 else
515 {
518 ->second;
519 }
520 }
521
522 /**
523 * @brief determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each
524 * quad point for each cell.
525 */
526 const auto &
528 {
529 if constexpr (std::is_same<ValueTypeBasisCoeff,
530 ValueTypeBasisData>::value)
531 {
532 return d_JxWData.find(d_quadratureID)->second;
533 }
534 else
535 {
536 return d_JxWBasisData.find(d_quadratureID)->second;
537 }
538 }
539
540 /**
541 * @brief Cell level stiffness matrix in ValueTypeBasisCoeff
542 */
543 const auto &
545 {
546 if constexpr (std::is_same<ValueTypeBasisCoeff,
547 ValueTypeBasisData>::value)
548 {
550 }
551 else
552 {
554 }
555 }
556
557
558 /**
559 * @brief Cell level stiffness matrix in ValueTypeBasisData
560 */
563
564
565 /**
566 * @brief Cell level mass matrix in ValueTypeBasisCoeff
567 */
568 const auto &
570 {
571 if constexpr (std::is_same<ValueTypeBasisCoeff,
572 ValueTypeBasisData>::value)
573 {
575 }
576 else
577 {
579 }
580 }
581
582
583 /**
584 * @brief Cell level mass matrix in ValueTypeBasisData
585 */
588
589
590 /**
591 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
592 */
593 const auto &
595 {
596 if constexpr (std::is_same<ValueTypeBasisCoeff,
597 ValueTypeBasisData>::value)
598 {
600 }
601 else
602 {
604 }
605 }
606
607 /**
608 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff
609 */
610 const auto &
612 {
613 if constexpr (std::is_same<ValueTypeBasisCoeff,
614 ValueTypeBasisData>::value)
615 {
617 }
618 else
619 {
621 }
622 }
623
624 /**
625 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff
626 */
627 const auto &
629 {
630 if constexpr (std::is_same<ValueTypeBasisCoeff,
631 ValueTypeBasisData>::value)
632 {
634 }
635 else
636 {
638 }
639 }
640
641 /**
642 * @brief Cell level diagonal mass matrix in ValueTypeBasisCoeff
643 */
644 const auto &
646 {
647 if constexpr (std::is_same<ValueTypeBasisCoeff,
648 ValueTypeBasisData>::value)
649 {
651 }
652 else
653 {
655 }
656 }
657
658
659 /**
660 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
661 */
662 const auto &
664 {
665 if constexpr (std::is_same<ValueTypeBasisCoeff,
666 ValueTypeBasisData>::value)
667 {
669 }
670 else
671 {
673 }
674 }
675
676
677
678 /**
679 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
680 */
681 const auto &
683 {
684 if constexpr (std::is_same<ValueTypeBasisCoeff,
685 ValueTypeBasisData>::value)
686 {
688 }
689 else
690 {
692 }
693 }
694
695
696 /**
697 * @brief sqrt diagonal mass matrix in ValueTypeBasisCoeff
698 */
699 const auto &
701 {
702 if constexpr (std::is_same<ValueTypeBasisCoeff,
703 ValueTypeBasisData>::value)
704 {
706 }
707 else
708 {
710 }
711 }
712
713 /**
714 * @brief diagonal mass matrix in ValueTypeBasisCoeff
715 */
716 const auto &
718 {
719 if constexpr (std::is_same<ValueTypeBasisCoeff,
720 ValueTypeBasisData>::value)
721 {
723 }
724 else
725 {
727 }
728 }
729
730
731 /**
732 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData
733 */
736
737
738 /**
739 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
740 */
743
744 /**
745 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
746 */
749 memorySpace> &
751
752 /**
753 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisData
754 */
757
758
759 /**
760 * @brief Cell level diagonal mass matrix in ValueTypeBasisData
761 */
764
765 /**
766 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisData
767 */
770
771 /**
772 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
773 */
776
777 /**
778 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
779 */
782 memorySpace> &
784
785 /**
786 * @brief sqrt diagonal mass matrix in ValueTypeBasisData
787 */
790
791 /**
792 * @brief diagonal mass matrix in ValueTypeBasisData
793 */
796
797 /**
798 * @brief diagonal stiffness matrix in ValueTypeBasisData
799 */
802
803 /**
804 * @brief diagonal inverse stiffness matrix in ValueTypeBasisData
805 */
808
809 /**
810 * @brief diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
811 */
814
815 /**
816 * @brief diagonal sqrt stiffness matrix in ValueTypeBasisData
817 */
820
821 /**
822 * @brief returns 2 if all cells on current processor are Cartesian,
823 * 1 if all cells on current processor are affine and 0 otherwise.
824 */
827
828 /**
829 * @brief returns the deal.ii cellID corresponing to given cell Index.
830 * @param[in] iElem cell Index
831 */
832 dealii::CellId
833 cellID(const dftfe::uInt iElem) const;
834 /**
835 * @brief returns the deal.ii cell_iterator corresponing to given cell Index.
836 * @param[in] iElem cell Index
837 */
838
839 dealii::DoFHandler<3>::active_cell_iterator
840 getCellIterator(const dftfe::uInt iElem) const;
841
842 /**
843 * @brief returns the cell index corresponding to given deal.ii cellID.
844 * @param[in] iElem cell Index
845 */
847 cellIndex(const dealii::CellId cellid) const;
848
849 /**
850 * @brief Creates a multivector.
851 * @param[in] blocksize Number of vectors in the multivector.
852 * @param[out] multiVector the created multivector.
853 */
854 void
856 const dftfe::uInt blocksize,
858 &multiVector) const;
859
860 /**
861 * @brief Creates a multivector.
862 * @param[in] blocksize Number of vectors in the multivector.
863 * @param[out] multiVector the created multivector.
864 */
865 void
867 const dftfe::uInt blocksize,
870 memorySpace> &multiVector) const;
871
872 /**
873 * @brief Creates scratch multivectors.
874 * @param[in] vecBlockSize Number of vectors in the multivector.
875 * @param[out] numMultiVecs number of scratch multivectors needed with
876 * this vecBlockSize.
877 */
878 void
880 const dftfe::uInt numMultiVecs = 1) const;
881
882 /**
883 * @brief Creates single precision scratch multivectors.
884 * @param[in] vecBlockSize Number of vectors in the multivector.
885 * @param[out] numMultiVecs number of scratch multivectors needed with
886 * this vecBlockSize.
887 */
888 void
890 const dftfe::uInt vecBlockSize,
891 const dftfe::uInt numMultiVecs = 1) const;
892
893 /**
894 * @brief Clears scratch multivectors.
895 */
896 void
898
899 /**
900 * @brief Gets scratch multivectors.
901 * @param[in] vecBlockSize Number of vectors in the multivector.
902 * @param[out] numMultiVecs index of the multivector among those with the
903 * same vecBlockSize.
904 */
906 getMultiVector(const dftfe::uInt vecBlockSize,
907 const dftfe::uInt index = 0) const;
908
909 /**
910 * @brief Gets single precision scratch multivectors.
911 * @param[in] vecBlockSize Number of vectors in the multivector.
912 * @param[out] numMultiVecs index of the multivector among those with the
913 * same vecBlockSize.
914 */
917 memorySpace> &
919 const dftfe::uInt index = 0) const;
920
921 /**
922 * @brief Apply constraints on given multivector.
923 * @param[inout] multiVector the given multivector.
924 */
925 void
927 memorySpace> &multiVector,
928 dftfe::uInt constraintIndex =
929 std::numeric_limits<dftfe::uInt>::max()) const;
930
931
932
933 /**
934 * @brief Return the underlying deal.II matrixfree object.
935 */
936 const dealii::MatrixFree<3, ValueTypeBasisData> &
938
939 /**
940 * @brief Return the underlying deal.II dofhandler object.
941 */
942 const dealii::DoFHandler<3> &
944
945
946
947 std::vector<dftUtils::constraintMatrixInfo<memorySpace>> d_constraintInfo;
949 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
951 const dealii::MatrixFree<3, ValueTypeBasisData> *d_matrixFreeDataPtr;
954 std::map<dftfe::uInt,
955 dftfe::utils::MemoryStorage<ValueTypeBasisData,
958 dftfe::utils::MemoryStorage<ValueTypeBasisData,
963 std::vector<dealii::CellId> d_cellIndexToCellIdMap;
964 std::vector<dealii::DoFHandler<3>::active_cell_iterator>
966 std::map<dealii::CellId, dftfe::uInt> d_cellIdToCellIndexMap;
967 std::map<dftfe::uInt,
970 std::map<dftfe::uInt,
973 std::map<dftfe::uInt,
976 std::map<dftfe::uInt,
979 std::map<dftfe::uInt,
982 std::map<dftfe::uInt,
985 std::map<dftfe::uInt,
988 std::map<dftfe::uInt,
991
992 std::map<dftfe::uInt,
995 std::map<dftfe::uInt,
998 std::map<dftfe::uInt,
1001 std::map<dftfe::uInt,
1004 std::map<dftfe::uInt,
1007 std::map<dftfe::uInt,
1010 std::map<dftfe::uInt,
1013
1024
1027 memorySpace>
1051 memorySpace>
1063 mutable std::map<
1065 std::vector<
1068
1085
1102 std::map<std::pair<dftfe::uInt, dftfe::uInt>,
1105 mutable std::map<
1109 memorySpace>>>
1111
1112 std::vector<dftfe::uInt> d_quadratureIDsVector;
1115 std::vector<dftfe::uInt> d_nQuadsPerCell;
1125 std::vector<UpdateFlags> d_updateFlags;
1126
1127 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
1129
1130
1131 void
1133 distributedCPUVec<double> &nodalField,
1134 const dftfe::uInt dofHandlerId,
1135 const dftfe::uInt quadratureId,
1137 &quadratureValueData,
1139 &quadratureGradValueData,
1141 &quadratureHessianValueData,
1142 const bool isEvaluateGradData = false,
1143 const bool isEvaluateHessianData = false,
1144 const bool isEvaluateData = true) const;
1145
1146 void
1148 const distributedCPUVec<double> &nodalField,
1149 const dftfe::uInt dofHandlerId,
1150 const dftfe::uInt quadratureId,
1152 &quadratureValueData,
1154 &quadratureGradValueData,
1156 &quadratureHessianValueData,
1157 const bool isEvaluateGradData = false,
1158 const bool isEvaluateHessianData = false,
1159 const bool isEvaluateData = true) const;
1160
1161 /**
1162 * @brief Interpolate process level nodal data to cell level quadrature data.
1163 * @param[in] nodalData process level nodal data, the multivector should
1164 * already have ghost data and constraints should have been applied.
1165 * @param[out] quadratureValues Cell level quadrature values, indexed by
1166 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1167 * @param[out] quadratureGradients Cell level quadrature gradients,
1168 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1169 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1170 */
1171 void
1173 memorySpace> &nodalData,
1174 ValueTypeBasisCoeff *quadratureValues,
1175 ValueTypeBasisCoeff *quadratureGradients = NULL) const;
1176
1177 // FIXME Untested function
1178 /**
1179 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1180 * @param[in] quadratureValues Cell level quadrature values, indexed by
1181 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1182 * @param[in] quadratureGradients Cell level quadrature gradients,
1183 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1184 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1185 * @param[out] nodalData process level nodal data.
1186 */
1187 void
1189 ValueTypeBasisCoeff *quadratureValues,
1190 ValueTypeBasisCoeff *quadratureGradients,
1192 &nodalData,
1194 &mapQuadIdToProcId) const;
1195
1196 /**
1197 * @brief Get cell level nodal data from process level nodal data.
1198 * @param[in] nodalData process level nodal data, the multivector should
1199 * already have ghost data and constraints should have been applied.
1200 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1201 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1202 */
1203 void
1206 &nodalData,
1207 ValueTypeBasisCoeff *cellNodalDataPtr) const;
1208 // FIXME Untested function
1209 /**
1210 * @brief Accumulate cell level nodal data into process level nodal data.
1211 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1212 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1213 * @param[out] nodalData process level nodal data.
1214 */
1215 void
1217 const ValueTypeBasisCoeff *cellNodalDataPtr,
1219 &nodalData) const;
1220
1221 /**
1222 * @brief Interpolate process level nodal data to cell level quadrature data.
1223 * @param[in] nodalData process level nodal data, the multivector should
1224 * already have ghost data and constraints should have been applied.
1225 * @param[out] quadratureValues Cell level quadrature values, indexed by
1226 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1227 * @param[out] quadratureGradients Cell level quadrature gradients,
1228 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1229 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1230 * @param[in] cellRange the range of cells for which interpolation has to
1231 * be done.
1232 */
1233 void
1235 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1236 memorySpace> &nodalData,
1237 ValueTypeBasisCoeff *quadratureValues,
1238 ValueTypeBasisCoeff *quadratureGradients,
1239 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1240
1241 /**
1242 * @brief Interpolate cell level nodal data to cell level quadrature data.
1243 * @param[in] nodalData cell level nodal data, the multivector should
1244 * already have ghost data and constraints should have been applied.
1245 * @param[out] quadratureValues Cell level quadrature values, indexed by
1246 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1247 * @param[out] quadratureGradients Cell level quadrature gradients,
1248 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1249 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1250 * @param[in] cellRange the range of cells for which interpolation has to
1251 * be done.
1252 */
1253 void
1255 const ValueTypeBasisCoeff *nodalData,
1256 ValueTypeBasisCoeff *quadratureValues,
1257 ValueTypeBasisCoeff *quadratureGradients,
1258 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1259
1260
1261 void
1262 interpolateQ1ToQ2(const double *Q1Field,
1263 const dftfe::uInt quadRule1,
1264 const dftfe::uInt quadRule2,
1265 double *Q2Field,
1266 const dftfe::uInt numComponents) const;
1267
1268 // FIXME Untested function
1269 /**
1270 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1271 * @param[in] quadratureValues Cell level quadrature values, indexed by
1272 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1273 * @param[in] quadratureGradients Cell level quadrature gradients,
1274 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1275 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1276 * @param[out] nodalData process level nodal data.
1277 * @param[in] cellRange the range of cells for which integration has to be
1278 * done.
1279 */
1280 void
1282 const ValueTypeBasisCoeff *quadratureValues,
1283 const ValueTypeBasisCoeff *quadratureGradients,
1285 &nodalData,
1287 &mapQuadIdToProcId,
1288 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1289
1290
1291 /**
1292 * @brief Get cell level nodal data from process level nodal data.
1293 * @param[in] nodalData process level nodal data, the multivector should
1294 * already have ghost data and constraints should have been applied.
1295 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1296 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1297 * @param[in] cellRange the range of cells for which extraction has to be
1298 * done.
1299 */
1300 void
1302 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1303 memorySpace> &nodalData,
1304 ValueTypeBasisCoeff *cellNodalDataPtr,
1305 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1306
1307 // FIXME Untested function
1308 /**
1309 * @brief Accumulate cell level nodal data into process level nodal data.
1310 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1311 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1312 * @param[out] nodalData process level nodal data.
1313 * @param[in] cellRange the range of cells for which extraction has to be
1314 * done.
1315 */
1316 void
1318 const ValueTypeBasisCoeff *cellNodalDataPtr,
1320 &nodalData,
1321 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1322 };
1323 } // end of namespace basis
1324} // end of namespace dftfe
1325
1326#endif // dftfeBasisOperations_h
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