DFT-FE 1.1.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 };
47
48 inline UpdateFlags
49 operator|(const UpdateFlags f1, const UpdateFlags f2)
50 {
51 return static_cast<UpdateFlags>(static_cast<dftfe::uInt>(f1) |
52 static_cast<dftfe::uInt>(f2));
53 }
54
55
56
57 inline UpdateFlags &
59 {
60 f1 = f1 | f2;
61 return f1;
62 }
63
64
65 inline UpdateFlags
66 operator&(const UpdateFlags f1, const UpdateFlags f2)
67 {
68 return static_cast<UpdateFlags>(static_cast<dftfe::uInt>(f1) &
69 static_cast<dftfe::uInt>(f2));
70 }
71
72
73 inline UpdateFlags &
75 {
76 f1 = f1 & f2;
77 return f1;
78 }
79
80
81 template <typename ValueTypeBasisCoeff,
82 typename ValueTypeBasisData,
83 dftfe::utils::MemorySpace memorySpace>
85 {
86 protected:
97 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
99
100 public:
101 /**
102 * @brief Constructor
103 */
106 BLASWrapperPtr);
107
108
109 /**
110 * @brief Default Destructor
111 */
113
114 /**
115 * @brief Clears the FEBasisOperations internal storage.
116 */
117 void
119
120 /**
121 * @brief fills required data structures for the given dofHandlerID
122 * @param[in] matrixFreeData MatrixFree object.
123 * @param[in] constraintsVector std::vector of AffineConstraints, should
124 * be the same vector which was passed for the construction of the given
125 * MatrixFree object.
126 * @param[in] dofHandlerID dofHandler index to be used for getting data
127 * from the MatrixFree object.
128 * @param[in] quadratureID std::vector of quadratureIDs to be used, should
129 * be the same IDs which were used during the construction of the given
130 * MatrixFree object.
131 */
132 void
133 init(dealii::MatrixFree<3, ValueTypeBasisData> &matrixFreeData,
134 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
135 &constraintsVector,
136 const dftfe::uInt &dofHandlerID,
137 const std::vector<dftfe::uInt> &quadratureID,
138 const std::vector<UpdateFlags> updateFlags);
139
140 /**
141 * @brief fills required data structures from another FEBasisOperations object
142 * @param[in] basisOperationsSrc Source FEBasisOperations object.
143 */
144 template <dftfe::utils::MemorySpace memorySpaceSrc>
145 void
146 init(const FEBasisOperations<ValueTypeBasisCoeff,
147 ValueTypeBasisData,
148 memorySpaceSrc> &basisOperationsSrc);
149
150 /**
151 * @brief sets internal variables and optionally resizes internal temp storage for interpolation operations
152 * @param[in] vecBlockSize block size to used for operations on vectors,
153 * this has to be set to the exact value before any such operations are
154 * called.
155 * @param[in] cellBlockSize block size to used for cells, this has to be
156 * set to a value greater than or equal to the required value before any
157 * such operations are called
158 * @param[in] quadratureID Quadrature index to be used.
159 * @param[in] isResizeTempStorage whether to resize internal tempstorage.
160 */
161 void
162 reinit(const dftfe::uInt &vecBlockSize,
163 const dftfe::uInt &cellBlockSize,
164 const dftfe::uInt &quadratureID,
165 const bool isResizeTempStorageForInerpolation = true,
166 const bool isResizeTempStorageForCellMatrices = false);
167
171
172 // private:
173
174
175 /**
176 * @brief Initializes indexset maps from process level indices to cell level indices for a single vector, also initializes cell index to cellid map.
177 */
178 void
180
181 /**
182 * @brief Initializes indexset maps from process level indices to cell level indices for multivectors.
183 */
184 void
186
187 /**
188 * @brief Initializes the constraintMatrixInfo object.
189 */
190 void
192
193 /**
194 * @brief Reinitializes the constraintMatrixInfo object.
195 */
196 void
198 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
199 &constraintsVector);
200
201 /**
202 * @brief Constructs the MPIPatternP2P object.
203 */
204 void
206
207 /**
208 * @brief Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
209 */
210 void
212
213 /**
214 * @brief Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
215 */
216 void
218
219
220 /**
221 * @brief Computes the cell-level stiffness matrix.
222 */
223 void
225 const dftfe::uInt cellsBlockSize,
226 const bool basisType = false,
227 const bool ceoffType = true);
228
229 void
231 const dftfe::uInt cellsBlockSize,
232 const bool basisType = false,
233 const bool ceoffType = true);
234
235 void
237 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
240 &weightedCellMassMatrix) const;
241
242 void
244 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
247 &weightedCellNjGradNiMatrix) const;
248
249 void
251 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
254 &weightedCellNjGradNiPlusNiGradNjMatrix) const;
255
256 void
258 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
261 &weightedCellNjGradNiPlusNiGradNjMatrix) const;
262
263 void
265 const std::pair<dftfe::uInt, dftfe::uInt> cellRangeTotal,
268 &weightedCellStiffnessMatrix) const;
269
270 void
271 computeInverseSqrtMassVector(const bool basisType = true,
272 const bool ceoffType = false);
273
274 /**
275 * @brief Computes the stiffness matrix
276 * \grad Ni . \grad Ni.
277 */
278 void
279 computeStiffnessVector(const bool basisType = true,
280 const bool ceoffType = false);
281
282 /**
283 * @brief Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in reinit.
284 */
285 void
286 resizeTempStorage(const bool isResizeTempStorageForInerpolation,
287 const bool isResizeTempStorageForCellMatrices);
288
289 /**
290 * @brief Number of quadrature points per cell for the quadratureID set in reinit.
291 */
294
295 /**
296 * @brief Number of vectors set in reinit.
297 */
299 nVectors() const;
300 /**
301 * @brief Number of DoFs per cell for the dofHandlerID set in init.
302 */
305
306 /**
307 * @brief Number of locally owned cells on the current processor.
308 */
310 nCells() const;
311
312 /**
313 * @brief Number of DoFs on the current processor, locally owned + ghosts.
314 */
317
318 /**
319 * @brief Number of locally owned DoFs on the current processor.
320 */
322 nOwnedDofs() const;
323
324 /**
325 * @brief Shape function values at quadrature points.
326 * @param[in] transpose if false the the data is indexed as [iQuad *
327 * d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
328 * d_nQuadsPerCell + iQuad].
329 */
331 shapeFunctionData(bool transpose = false) const;
332
333 /**
334 * @brief Shape function gradient values at quadrature points.
335 * @param[in] transpose if false the the data is indexed as [iDim *
336 * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
337 * if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
338 * iNode * d_nQuadsPerCell + iQuad].
339 */
341 shapeFunctionGradientData(bool transpose = false) const;
342
343 /**
344 * @brief Inverse Jacobian matrices, for cartesian cells returns the
345 * diagonal elements of the inverse Jacobian matrices for each cell, for
346 * affine cells returns the 3x3 inverse Jacobians for each cell otherwise
347 * returns the 3x3 inverse Jacobians at each quad point for each cell.
348 */
351
352 /**
353 * @brief determinant of Jacobian times the quadrature weight at each
354 * quad point for each cell.
355 */
357 JxW() const;
358
359 /**
360 * @brief quad point coordinates for each cell.
361 */
362 const dftfe::utils::MemoryStorage<ValueTypeBasisData,
364 quadPoints() const;
365
366 /**
367 * @brief Shape function values at quadrature points in ValueTypeBasisData.
368 * @param[in] transpose if false the the data is indexed as [iQuad *
369 * d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
370 * d_nQuadsPerCell + iQuad].
371 */
372 const auto &
373 shapeFunctionBasisData(bool transpose = false) const
374 {
375 if constexpr (std::is_same<ValueTypeBasisCoeff,
376 ValueTypeBasisData>::value)
377 {
378 return transpose ?
381 }
382 else
383 {
384 return transpose ?
386 ->second :
388 }
389 }
390 /**
391 * @brief Shape function gradient values at quadrature points in ValueTypeBasisData.
392 * @param[in] transpose if false the the data is indexed as [iDim *
393 * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
394 * if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
395 * iNode * d_nQuadsPerCell + iQuad].
396 */
397 const auto &
398 shapeFunctionGradientBasisData(bool transpose = false) const
399 {
400 if constexpr (std::is_same<ValueTypeBasisCoeff,
401 ValueTypeBasisData>::value)
402 {
403 return transpose ?
405 ->second :
407 }
408 else
409 {
410 return transpose ?
412 .find(d_quadratureID)
413 ->second :
415 ->second;
416 }
417 }
418
419 /**
420 * @brief Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the
421 * diagonal elements of the inverse Jacobian matrices for each cell, for
422 * affine cells returns the 3x3 inverse Jacobians for each cell otherwise
423 * returns the 3x3 inverse Jacobians at each quad point for each cell.
424 */
425 const auto &
427 {
428 if constexpr (std::is_same<ValueTypeBasisCoeff,
429 ValueTypeBasisData>::value)
430 {
433 ->second;
434 }
435 else
436 {
439 ->second;
440 }
441 }
442
443 /**
444 * @brief determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each
445 * quad point for each cell.
446 */
447 const auto &
449 {
450 if constexpr (std::is_same<ValueTypeBasisCoeff,
451 ValueTypeBasisData>::value)
452 {
453 return d_JxWData.find(d_quadratureID)->second;
454 }
455 else
456 {
457 return d_JxWBasisData.find(d_quadratureID)->second;
458 }
459 }
460
461 /**
462 * @brief Cell level stiffness matrix in ValueTypeBasisCoeff
463 */
464 const auto &
466 {
467 if constexpr (std::is_same<ValueTypeBasisCoeff,
468 ValueTypeBasisData>::value)
469 {
471 }
472 else
473 {
475 }
476 }
477
478
479 /**
480 * @brief Cell level stiffness matrix in ValueTypeBasisData
481 */
484
485
486 /**
487 * @brief Cell level mass matrix in ValueTypeBasisCoeff
488 */
489 const auto &
491 {
492 if constexpr (std::is_same<ValueTypeBasisCoeff,
493 ValueTypeBasisData>::value)
494 {
496 }
497 else
498 {
500 }
501 }
502
503
504 /**
505 * @brief Cell level mass matrix in ValueTypeBasisData
506 */
509
510
511 /**
512 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
513 */
514 const auto &
516 {
517 if constexpr (std::is_same<ValueTypeBasisCoeff,
518 ValueTypeBasisData>::value)
519 {
521 }
522 else
523 {
525 }
526 }
527
528 /**
529 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff
530 */
531 const auto &
533 {
534 if constexpr (std::is_same<ValueTypeBasisCoeff,
535 ValueTypeBasisData>::value)
536 {
538 }
539 else
540 {
542 }
543 }
544
545 /**
546 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff
547 */
548 const auto &
550 {
551 if constexpr (std::is_same<ValueTypeBasisCoeff,
552 ValueTypeBasisData>::value)
553 {
555 }
556 else
557 {
559 }
560 }
561
562 /**
563 * @brief Cell level diagonal mass matrix in ValueTypeBasisCoeff
564 */
565 const auto &
567 {
568 if constexpr (std::is_same<ValueTypeBasisCoeff,
569 ValueTypeBasisData>::value)
570 {
572 }
573 else
574 {
576 }
577 }
578
579
580 /**
581 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
582 */
583 const auto &
585 {
586 if constexpr (std::is_same<ValueTypeBasisCoeff,
587 ValueTypeBasisData>::value)
588 {
590 }
591 else
592 {
594 }
595 }
596
597
598
599 /**
600 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
601 */
602 const auto &
604 {
605 if constexpr (std::is_same<ValueTypeBasisCoeff,
606 ValueTypeBasisData>::value)
607 {
609 }
610 else
611 {
613 }
614 }
615
616
617 /**
618 * @brief sqrt diagonal mass matrix in ValueTypeBasisCoeff
619 */
620 const auto &
622 {
623 if constexpr (std::is_same<ValueTypeBasisCoeff,
624 ValueTypeBasisData>::value)
625 {
627 }
628 else
629 {
631 }
632 }
633
634 /**
635 * @brief diagonal mass matrix in ValueTypeBasisCoeff
636 */
637 const auto &
639 {
640 if constexpr (std::is_same<ValueTypeBasisCoeff,
641 ValueTypeBasisData>::value)
642 {
644 }
645 else
646 {
648 }
649 }
650
651
652 /**
653 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData
654 */
657
658
659 /**
660 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
661 */
664
665 /**
666 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
667 */
670 memorySpace> &
672
673 /**
674 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisData
675 */
678
679
680 /**
681 * @brief Cell level diagonal mass matrix in ValueTypeBasisData
682 */
685
686 /**
687 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisData
688 */
691
692 /**
693 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
694 */
697
698 /**
699 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
700 */
703 memorySpace> &
705
706 /**
707 * @brief sqrt diagonal mass matrix in ValueTypeBasisData
708 */
711
712 /**
713 * @brief diagonal mass matrix in ValueTypeBasisData
714 */
717
718 /**
719 * @brief diagonal stiffness matrix in ValueTypeBasisData
720 */
723
724 /**
725 * @brief diagonal inverse stiffness matrix in ValueTypeBasisData
726 */
729
730 /**
731 * @brief diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
732 */
735
736 /**
737 * @brief diagonal sqrt stiffness matrix in ValueTypeBasisData
738 */
741
742 /**
743 * @brief returns 2 if all cells on current processor are Cartesian,
744 * 1 if all cells on current processor are affine and 0 otherwise.
745 */
748
749 /**
750 * @brief returns the deal.ii cellID corresponing to given cell Index.
751 * @param[in] iElem cell Index
752 */
753 dealii::CellId
754 cellID(const dftfe::uInt iElem) const;
755 /**
756 * @brief returns the deal.ii cell_iterator corresponing to given cell Index.
757 * @param[in] iElem cell Index
758 */
759
760 dealii::DoFHandler<3>::active_cell_iterator
761 getCellIterator(const dftfe::uInt iElem) const;
762
763 /**
764 * @brief returns the cell index corresponding to given deal.ii cellID.
765 * @param[in] iElem cell Index
766 */
768 cellIndex(const dealii::CellId cellid) const;
769
770 /**
771 * @brief Creates a multivector.
772 * @param[in] blocksize Number of vectors in the multivector.
773 * @param[out] multiVector the created multivector.
774 */
775 void
777 const dftfe::uInt blocksize,
779 &multiVector) const;
780
781 /**
782 * @brief Creates a multivector.
783 * @param[in] blocksize Number of vectors in the multivector.
784 * @param[out] multiVector the created multivector.
785 */
786 void
788 const dftfe::uInt blocksize,
791 memorySpace> &multiVector) const;
792
793 /**
794 * @brief Creates scratch multivectors.
795 * @param[in] vecBlockSize Number of vectors in the multivector.
796 * @param[out] numMultiVecs number of scratch multivectors needed with
797 * this vecBlockSize.
798 */
799 void
801 const dftfe::uInt numMultiVecs = 1) const;
802
803 /**
804 * @brief Creates single precision scratch multivectors.
805 * @param[in] vecBlockSize Number of vectors in the multivector.
806 * @param[out] numMultiVecs number of scratch multivectors needed with
807 * this vecBlockSize.
808 */
809 void
811 const dftfe::uInt vecBlockSize,
812 const dftfe::uInt numMultiVecs = 1) const;
813
814 /**
815 * @brief Clears scratch multivectors.
816 */
817 void
819
820 /**
821 * @brief Gets scratch multivectors.
822 * @param[in] vecBlockSize Number of vectors in the multivector.
823 * @param[out] numMultiVecs index of the multivector among those with the
824 * same vecBlockSize.
825 */
827 getMultiVector(const dftfe::uInt vecBlockSize,
828 const dftfe::uInt index = 0) const;
829
830 /**
831 * @brief Gets single precision scratch multivectors.
832 * @param[in] vecBlockSize Number of vectors in the multivector.
833 * @param[out] numMultiVecs index of the multivector among those with the
834 * same vecBlockSize.
835 */
838 memorySpace> &
840 const dftfe::uInt index = 0) const;
841
842 /**
843 * @brief Apply constraints on given multivector.
844 * @param[inout] multiVector the given multivector.
845 */
846 void
848 memorySpace> &multiVector,
849 dftfe::uInt constraintIndex =
850 std::numeric_limits<dftfe::uInt>::max()) const;
851
852
853
854 /**
855 * @brief Return the underlying deal.II matrixfree object.
856 */
857 const dealii::MatrixFree<3, ValueTypeBasisData> &
859
860 /**
861 * @brief Return the underlying deal.II dofhandler object.
862 */
863 const dealii::DoFHandler<3> &
865
866
867
868 std::vector<dftUtils::constraintMatrixInfo<memorySpace>> d_constraintInfo;
870 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
872 const dealii::MatrixFree<3, ValueTypeBasisData> *d_matrixFreeDataPtr;
875 std::map<dftfe::uInt,
876 dftfe::utils::MemoryStorage<ValueTypeBasisData,
881 std::vector<dealii::CellId> d_cellIndexToCellIdMap;
882 std::vector<dealii::DoFHandler<3>::active_cell_iterator>
884 std::map<dealii::CellId, dftfe::uInt> d_cellIdToCellIndexMap;
885 std::map<dftfe::uInt,
888 std::map<dftfe::uInt,
891 std::map<dftfe::uInt,
894 std::map<dftfe::uInt,
897 std::map<dftfe::uInt,
900 std::map<dftfe::uInt,
903 std::map<dftfe::uInt,
906
907 std::map<dftfe::uInt,
910 std::map<dftfe::uInt,
913 std::map<dftfe::uInt,
916 std::map<dftfe::uInt,
919 std::map<dftfe::uInt,
922 std::map<dftfe::uInt,
925
936
939 memorySpace>
963 memorySpace>
975 mutable std::map<
977 std::vector<
980
997
1014
1015
1016 mutable std::map<
1020 memorySpace>>>
1022
1023 std::vector<dftfe::uInt> d_quadratureIDsVector;
1026 std::vector<dftfe::uInt> d_nQuadsPerCell;
1036 std::vector<UpdateFlags> d_updateFlags;
1037
1038 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
1040
1041
1042 /**
1043 * @brief Interpolate process level nodal data to cell level quadrature data.
1044 * @param[in] nodalData process level nodal data, the multivector should
1045 * already have ghost data and constraints should have been applied.
1046 * @param[out] quadratureValues Cell level quadrature values, indexed by
1047 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1048 * @param[out] quadratureGradients Cell level quadrature gradients,
1049 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1050 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1051 */
1052 void
1054 memorySpace> &nodalData,
1055 ValueTypeBasisCoeff *quadratureValues,
1056 ValueTypeBasisCoeff *quadratureGradients = NULL) const;
1057
1058 // FIXME Untested function
1059 /**
1060 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1061 * @param[in] quadratureValues Cell level quadrature values, indexed by
1062 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1063 * @param[in] quadratureGradients Cell level quadrature gradients,
1064 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1065 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1066 * @param[out] nodalData process level nodal data.
1067 */
1068 void
1070 ValueTypeBasisCoeff *quadratureValues,
1071 ValueTypeBasisCoeff *quadratureGradients,
1073 &nodalData,
1075 &mapQuadIdToProcId) const;
1076
1077 /**
1078 * @brief Get cell level nodal data from process level nodal data.
1079 * @param[in] nodalData process level nodal data, the multivector should
1080 * already have ghost data and constraints should have been applied.
1081 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1082 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1083 */
1084 void
1087 &nodalData,
1088 ValueTypeBasisCoeff *cellNodalDataPtr) const;
1089 // FIXME Untested function
1090 /**
1091 * @brief Accumulate cell level nodal data into process level nodal data.
1092 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1093 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1094 * @param[out] nodalData process level nodal data.
1095 */
1096 void
1098 const ValueTypeBasisCoeff *cellNodalDataPtr,
1100 &nodalData) const;
1101
1102 /**
1103 * @brief Interpolate process level nodal data to cell level quadrature data.
1104 * @param[in] nodalData process level nodal data, the multivector should
1105 * already have ghost data and constraints should have been applied.
1106 * @param[out] quadratureValues Cell level quadrature values, indexed by
1107 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1108 * @param[out] quadratureGradients Cell level quadrature gradients,
1109 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1110 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1111 * @param[in] cellRange the range of cells for which interpolation has to
1112 * be done.
1113 */
1114 void
1116 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1117 memorySpace> &nodalData,
1118 ValueTypeBasisCoeff *quadratureValues,
1119 ValueTypeBasisCoeff *quadratureGradients,
1120 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1121
1122 /**
1123 * @brief Interpolate cell level nodal data to cell level quadrature data.
1124 * @param[in] nodalData cell level nodal data, the multivector should
1125 * already have ghost data and constraints should have been applied.
1126 * @param[out] quadratureValues Cell level quadrature values, indexed by
1127 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1128 * @param[out] quadratureGradients Cell level quadrature gradients,
1129 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1130 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1131 * @param[in] cellRange the range of cells for which interpolation has to
1132 * be done.
1133 */
1134 void
1136 const ValueTypeBasisCoeff *nodalData,
1137 ValueTypeBasisCoeff *quadratureValues,
1138 ValueTypeBasisCoeff *quadratureGradients,
1139 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1140
1141 // FIXME Untested function
1142 /**
1143 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1144 * @param[in] quadratureValues Cell level quadrature values, indexed by
1145 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1146 * @param[in] quadratureGradients Cell level quadrature gradients,
1147 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1148 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1149 * @param[out] nodalData process level nodal data.
1150 * @param[in] cellRange the range of cells for which integration has to be
1151 * done.
1152 */
1153 void
1155 const ValueTypeBasisCoeff *quadratureValues,
1156 const ValueTypeBasisCoeff *quadratureGradients,
1158 &nodalData,
1160 &mapQuadIdToProcId,
1161 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1162
1163
1164 /**
1165 * @brief Get cell level nodal data from process level nodal data.
1166 * @param[in] nodalData process level nodal data, the multivector should
1167 * already have ghost data and constraints should have been applied.
1168 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1169 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1170 * @param[in] cellRange the range of cells for which extraction has to be
1171 * done.
1172 */
1173 void
1175 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1176 memorySpace> &nodalData,
1177 ValueTypeBasisCoeff *cellNodalDataPtr,
1178 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1179
1180 // FIXME Untested function
1181 /**
1182 * @brief Accumulate cell level nodal data into process level nodal data.
1183 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1184 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1185 * @param[out] nodalData process level nodal data.
1186 * @param[in] cellRange the range of cells for which extraction has to be
1187 * done.
1188 */
1189 void
1191 const ValueTypeBasisCoeff *cellNodalDataPtr,
1193 &nodalData,
1194 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1195 };
1196 } // end of namespace basis
1197} // end of namespace dftfe
1198
1199#endif // dftfeBasisOperations_h
dftfe::uInt d_quadratureIndex
Definition FEBasisOperations.h:1025
dftfe::uInt nVectors() const
Number of vectors set in reinit.
void computeWeightedCellStiffnessMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellStiffnessMatrix) const
std::vector< dftfe::uInt > d_quadratureIDsVector
Definition FEBasisOperations.h:1023
const auto & massVector() const
diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:638
dealii::CellId cellID(const dftfe::uInt iElem) const
returns the deal.ii cellID corresponing to given cell Index.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > & quadPoints() const
quad point coordinates for each cell.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseSqrtMassVectorBasisData() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisData.
const auto & cellMassVector() const
Cell level diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:566
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseMassVectorCoeffType
Definition FEBasisOperations.h:966
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:986
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisData
Definition FEBasisOperations.h:915
const auto & inverseSqrtMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:584
bool areAllCellsAffine
Definition FEBasisOperations.h:1034
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:970
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseStiffnessVectorBasisData() const
diagonal inverse stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1007
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_cellInverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:940
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtMassVectorCoeffType
Definition FEBasisOperations.h:954
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_inverseJacobianData
Definition FEBasisOperations.h:887
const dealii::MatrixFree< 3, ValueTypeBasisData > * d_matrixFreeDataPtr
Definition FEBasisOperations.h:872
const auto & cellInverseMassVector() const
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:532
void initializeMPIPattern()
Constructs the MPIPatternP2P object.
void interpolate(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients=NULL) const
Interpolate process level nodal data to cell level quadrature data.
std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > d_cellIndexToCellIteratorMap
Definition FEBasisOperations.h:883
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock
Definition FEBasisOperations.h:91
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseMassVectorBasisType
Definition FEBasisOperations.h:935
void accumulateFromCellNodalDataKernel(const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const
Accumulate cell level nodal data into process level nodal data.
dftfe::uInt d_locallyOwnedSize
Definition FEBasisOperations.h:1033
std::map< dealii::CellId, dftfe::uInt > d_cellIdToCellIndexMap
Definition FEBasisOperations.h:884
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > > d_quadPoints
Definition FEBasisOperations.h:878
const dealii::DoFHandler< 3 > & getDofHandler() const
Return the underlying deal.II dofhandler object.
void computeCellMassMatrix(const dftfe::uInt quadratureID, const dftfe::uInt cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1001
dftfe::uInt nOwnedDofs() const
Number of locally owned DoFs on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_massVectorCoeffType
Definition FEBasisOperations.h:958
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataInternalLayout
Definition FEBasisOperations.h:896
dealii::DoFHandler< 3 >::active_cell_iterator getCellIterator(const dftfe::uInt iElem) const
returns the deal.ii cell_iterator corresponing to given cell Index.
void computeWeightedCellNjGradNiMinusNiGradNjMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellStiffnessMatrixBasisData() const
Cell level stiffness matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellInverseMassVectorBasisData() const
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
void initializeConstraints()
Initializes the constraintMatrixInfo object.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataTranspose
Definition FEBasisOperations.h:905
~FEBasisOperations()=default
Default Destructor.
dftfe::uInt d_dofHandlerID
Definition FEBasisOperations.h:1027
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseStiffnessVectorBasisType
Definition FEBasisOperations.h:984
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionData
Definition FEBasisOperations.h:893
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:990
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:988
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & inverseJacobians() const
Inverse Jacobian matrices, for cartesian cells returns the diagonal elements of the inverse Jacobian ...
void extractToCellNodalData(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr) const
Get cell level nodal data from process level nodal data.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionDataTranspose
Definition FEBasisOperations.h:902
const auto & inverseJacobiansBasisData() const
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of...
Definition FEBasisOperations.h:426
const auto & cellStiffnessMatrix() const
Cell level stiffness matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:465
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
Definition FEBasisOperations.h:98
void resizeTempStorage(const bool isResizeTempStorageForInerpolation, const bool isResizeTempStorageForCellMatrices)
Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in re...
void computeCellStiffnessMatrix(const dftfe::uInt quadratureID, const dftfe::uInt cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
Computes the cell-level stiffness matrix.
dftfe::uInt d_nCells
Definition FEBasisOperations.h:1029
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > zeroIndexVec
Definition FEBasisOperations.h:94
const auto & shapeFunctionBasisData(bool transpose=false) const
Shape function values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:373
void reinit(const dftfe::uInt &vecBlockSize, const dftfe::uInt &cellBlockSize, const dftfe::uInt &quadratureID, const bool isResizeTempStorageForInerpolation=true, const bool isResizeTempStorageForCellMatrices=false)
sets internal variables and optionally resizes internal temp storage for interpolation operations
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_stiffnessVectorBasisType
Definition FEBasisOperations.h:996
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock2
Definition FEBasisOperations.h:91
void createMultiVector(const dftfe::uInt blocksize, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector) const
Creates a multivector.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:968
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > mpiPatternP2P
Definition FEBasisOperations.h:1039
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1009
void computeInverseSqrtMassVector(const bool basisType=true, const bool ceoffType=false)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:999
void createMultiVectorSinglePrec(const dftfe::uInt blocksize, dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > &multiVector) const
Creates a multivector.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1005
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseStiffnessVectorBasisType
Definition FEBasisOperations.h:994
dftfe::uInt d_cellsBlockSize
Definition FEBasisOperations.h:1030
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessVectorBasisType
Definition FEBasisOperations.h:982
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtMassVectorCoeffType
Definition FEBasisOperations.h:974
dftfe::uInt d_nVectors
Definition FEBasisOperations.h:1028
dftfe::uInt d_nOMPThreads
Definition FEBasisOperations.h:869
std::vector< dftfe::uInt > d_nQuadsPerCell
Definition FEBasisOperations.h:1026
dftfe::uInt nDofsPerCell() const
Number of DoFs per cell for the dofHandlerID set in init.
void accumulateFromCellNodalData(const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData) const
Accumulate cell level nodal data into process level nodal data.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & massVectorBasisData() const
diagonal mass matrix in ValueTypeBasisData
void reinitializeConstraints(std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector)
Reinitializes the constraintMatrixInfo object.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassMatrixCoeffType
Definition FEBasisOperations.h:933
void computeWeightedCellNjGradNiPlusNiGradNjMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtMassVectorBasisType
Definition FEBasisOperations.h:972
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & shapeFunctionGradientData(bool transpose=false) const
Shape function gradient values at quadrature points.
std::vector< dealii::CellId > d_cellIndexToCellIdMap
Definition FEBasisOperations.h:881
void clearScratchMultiVectors() const
Clears scratch multivectors.
std::vector< UpdateFlags > d_updateFlags
Definition FEBasisOperations.h:1036
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_inverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:964
dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > & getMultiVectorSinglePrec(const dftfe::uInt vecBlockSize, const dftfe::uInt index=0) const
Gets single precision scratch multivectors.
const auto & sqrtMassVector() const
sqrt diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:621
dftfe::uInt d_quadratureID
Definition FEBasisOperations.h:1024
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_JxWData
Definition FEBasisOperations.h:890
void createScratchMultiVectors(const dftfe::uInt vecBlockSize, const dftfe::uInt numMultiVecs=1) const
Creates scratch multivectors.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:1013
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & JxW() const
determinant of Jacobian times the quadrature weight at each quad point for each cell.
std::map< dftfe::uInt, std::vector< dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > > > scratchMultiVectorsSinglePrec
Definition FEBasisOperations.h:1021
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessMatrixBasisType
Definition FEBasisOperations.h:927
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > d_flattenedCellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:880
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & cellInverseMassVectorBasisDataSinglePrec() const
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessVectorCoeffType
Definition FEBasisOperations.h:1003
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:944
void clear()
Clears the FEBasisOperations internal storage.
void integrateWithBasisKernel(const ValueTypeBasisCoeff *quadratureValues, const ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > &mapQuadIdToProcId, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const
Integrate cell level quadrature data times shape functions to process level nodal data.
void interpolateKernel(const ValueTypeBasisCoeff *nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const
Interpolate cell level nodal data to cell level quadrature data.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientData
Definition FEBasisOperations.h:899
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellSqrtMassVectorBasisData() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisData.
void integrateWithBasis(ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > &mapQuadIdToProcId) const
Integrate cell level quadrature data times shape functions to process level nodal data.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseMassVectorBasisType
Definition FEBasisOperations.h:960
const auto & cellSqrtMassVector() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:549
void createScratchMultiVectorsSinglePrec(const dftfe::uInt vecBlockSize, const dftfe::uInt numMultiVecs=1) const
Creates single precision scratch multivectors.
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisDataTranspose
Definition FEBasisOperations.h:921
void computeStiffnessVector(const bool basisType=true, const bool ceoffType=false)
Computes the stiffness matrix \grad Ni . \grad Ni.
dftfe::uInt cellIndex(const dealii::CellId cellid) const
returns the cell index corresponding to given deal.ii cellID.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsData
Definition FEBasisOperations.h:88
std::map< dftfe::uInt, std::vector< dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > > > scratchMultiVectors
Definition FEBasisOperations.h:979
dftfe::uInt nCells() const
Number of locally owned cells on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellValuesBlockCoeff
Definition FEBasisOperations.h:96
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessMatrixCoeffType
Definition FEBasisOperations.h:929
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_stiffnessVectorCoeffType
Definition FEBasisOperations.h:1011
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellMassVectorBasisData() const
Cell level diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:946
std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > * d_constraintsVector
Definition FEBasisOperations.h:871
const auto & inverseMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:603
dftfe::uInt d_localSize
Definition FEBasisOperations.h:1032
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & inverseMassVectorBasisDataSinglePrec() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellValuesBlock
Definition FEBasisOperations.h:91
void computeWeightedCellMassMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellMassMatrix) const
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & getMultiVector(const dftfe::uInt vecBlockSize, const dftfe::uInt index=0) const
Gets scratch multivectors.
void initializeShapeFunctionAndJacobianData()
Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtMassVectorBasisType
Definition FEBasisOperations.h:952
void initializeShapeFunctionAndJacobianBasisData()
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseMassVectorCoeffType
Definition FEBasisOperations.h:942
dftfe::uInt nQuadsPerCell() const
Number of quadrature points per cell for the quadratureID set in reinit.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:992
dftfe::uInt cellsTypeFlag() const
returns 2 if all cells on current processor are Cartesian, 1 if all cells on current processor are af...
const auto & JxWBasisData() const
determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each quad point for each...
Definition FEBasisOperations.h:448
void interpolateKernel(const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const
Interpolate process level nodal data to cell level quadrature data.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtStiffnessVectorBasisData() const
diagonal sqrt stiffness matrix in ValueTypeBasisData
void init(const FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpaceSrc > &basisOperationsSrc)
fills required data structures from another FEBasisOperations object
void computeWeightedCellNjGradNiMatrix(const std::pair< dftfe::uInt, dftfe::uInt > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiMatrix) const
void init(dealii::MatrixFree< 3, ValueTypeBasisData > &matrixFreeData, std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector, const dftfe::uInt &dofHandlerID, const std::vector< dftfe::uInt > &quadratureID, const std::vector< UpdateFlags > updateFlags)
fills required data structures for the given dofHandlerID
bool areAllCellsCartesian
Definition FEBasisOperations.h:1035
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellNodalData
Definition FEBasisOperations.h:88
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassVectorBasisType
Definition FEBasisOperations.h:948
void initializeFlattenedIndexMaps()
Initializes indexset maps from process level indices to cell level indices for multivectors.
const auto & cellInverseSqrtMassVector() const
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:515
void initializeIndexMaps()
Initializes indexset maps from process level indices to cell level indices for a single vector,...
dftfe::utils::MemoryStorage< dftfe::uInt, dftfe::utils::MemorySpace::HOST > d_cellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:874
dftfe::uInt d_nDofsPerCell
Definition FEBasisOperations.h:1031
std::vector< dftUtils::constraintMatrixInfo< memorySpace > > d_constraintInfo
Definition FEBasisOperations.h:868
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseSqrtStiffnessVectorBasisData() const
diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisData
Definition FEBasisOperations.h:918
FEBasisOperations(std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr)
Constructor.
void extractToCellNodalDataKernel(const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr, const std::pair< dftfe::uInt, dftfe::uInt > cellRange) const
Get cell level nodal data from process level nodal data.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassVectorCoeffType
Definition FEBasisOperations.h:950
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_JxWBasisData
Definition FEBasisOperations.h:912
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseMassVectorBasisData() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellMatrixBlock
Definition FEBasisOperations.h:92
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtMassVectorBasisData() const
sqrt diagonal mass matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassMatrixBasisType
Definition FEBasisOperations.h:931
void distribute(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector, dftfe::uInt constraintIndex=std::numeric_limits< dftfe::uInt >::max()) const
Apply constraints on given multivector.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & stiffnessVectorBasisData() const
diagonal stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsDataNonAffine
Definition FEBasisOperations.h:89
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisDataTranspose
Definition FEBasisOperations.h:924
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & shapeFunctionData(bool transpose=false) const
Shape function values at quadrature points.
dftfe::uInt nRelaventDofs() const
Number of DoFs on the current processor, locally owned + ghosts.
const auto & cellMassMatrix() const
Cell level mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:490
const auto & shapeFunctionGradientBasisData(bool transpose=false) const
Shape function gradient values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:398
const dealii::MatrixFree< 3, ValueTypeBasisData > & matrixFreeData() const
Return the underlying deal.II matrixfree object.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellInverseSqrtMassVectorBasisData() const
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellMassMatrixBasisData() const
Cell level mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< dftfe::uInt, dftfe::utils::MemorySpace::HOST > & getFlattenedMapsHost()
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_massVectorBasisType
Definition FEBasisOperations.h:956
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_inverseJacobianBasisData
Definition FEBasisOperations.h:909
Definition BLASWrapper.h:35
An class template to encapsulate a MultiVector. A MultiVector is a collection of vectors belonging t...
Definition MultiVector.h:127
Definition MemoryStorage.h:33
Definition FEBasisOperations.h:30
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
Definition FEBasisOperations.h:49
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:58
UpdateFlags
Definition FEBasisOperations.h:32
@ update_values
Definition FEBasisOperations.h:35
@ update_quadpoints
Definition FEBasisOperations.h:41
@ update_transpose
Definition FEBasisOperations.h:39
@ update_inversejacobians
Definition FEBasisOperations.h:43
@ update_jxw
Definition FEBasisOperations.h:45
@ update_default
Definition FEBasisOperations.h:33
@ update_gradients
Definition FEBasisOperations.h:37
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
Definition FEBasisOperations.h:66
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:74
MemorySpace
Definition MemorySpaceType.h:33
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10
T type
Definition dftfeDataTypes.h:113