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<unsigned int>(f1) |
52 static_cast<unsigned int>(f2));
53 }
54
55
56
57 inline UpdateFlags &
59 {
60 f1 = f1 | f2;
61 return f1;
62 }
63
64
65 inline UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
66 {
67 return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) &
68 static_cast<unsigned int>(f2));
69 }
70
71
72 inline UpdateFlags &
74 {
75 f1 = f1 & f2;
76 return f1;
77 }
78
79
80 template <typename ValueTypeBasisCoeff,
81 typename ValueTypeBasisData,
82 dftfe::utils::MemorySpace memorySpace>
84 {
85 protected:
96 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
98
99 public:
100 /**
101 * @brief Constructor
102 */
105 BLASWrapperPtr);
106
107
108 /**
109 * @brief Default Destructor
110 */
112
113 /**
114 * @brief Clears the FEBasisOperations internal storage.
115 */
116 void
118
119 /**
120 * @brief fills required data structures for the given dofHandlerID
121 * @param[in] matrixFreeData MatrixFree object.
122 * @param[in] constraintsVector std::vector of AffineConstraints, should
123 * be the same vector which was passed for the construction of the given
124 * MatrixFree object.
125 * @param[in] dofHandlerID dofHandler index to be used for getting data
126 * from the MatrixFree object.
127 * @param[in] quadratureID std::vector of quadratureIDs to be used, should
128 * be the same IDs which were used during the construction of the given
129 * MatrixFree object.
130 */
131 void
132 init(dealii::MatrixFree<3, ValueTypeBasisData> &matrixFreeData,
133 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
134 & constraintsVector,
135 const unsigned int & dofHandlerID,
136 const std::vector<unsigned int> &quadratureID,
137 const std::vector<UpdateFlags> updateFlags);
138
139 /**
140 * @brief fills required data structures from another FEBasisOperations object
141 * @param[in] basisOperationsSrc Source FEBasisOperations object.
142 */
143 template <dftfe::utils::MemorySpace memorySpaceSrc>
144 void
145 init(const FEBasisOperations<ValueTypeBasisCoeff,
146 ValueTypeBasisData,
147 memorySpaceSrc> &basisOperationsSrc);
148
149 /**
150 * @brief sets internal variables and optionally resizes internal temp storage for interpolation operations
151 * @param[in] vecBlockSize block size to used for operations on vectors,
152 * this has to be set to the exact value before any such operations are
153 * called.
154 * @param[in] cellBlockSize block size to used for cells, this has to be
155 * set to a value greater than or equal to the required value before any
156 * such operations are called
157 * @param[in] quadratureID Quadrature index to be used.
158 * @param[in] isResizeTempStorage whether to resize internal tempstorage.
159 */
160 void
161 reinit(const unsigned int &vecBlockSize,
162 const unsigned int &cellBlockSize,
163 const unsigned int &quadratureID,
164 const bool isResizeTempStorageForInerpolation = true,
165 const bool isResizeTempStorageForCellMatrices = false);
166
170
171 // private:
172
173
174 /**
175 * @brief Initializes indexset maps from process level indices to cell level indices for a single vector, also initializes cell index to cellid map.
176 */
177 void
179
180 /**
181 * @brief Initializes indexset maps from process level indices to cell level indices for multivectors.
182 */
183 void
185
186 /**
187 * @brief Initializes the constraintMatrixInfo object.
188 */
189 void
191
192 /**
193 * @brief Reinitializes the constraintMatrixInfo object.
194 */
195 void
197 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
198 &constraintsVector);
199
200 /**
201 * @brief Constructs the MPIPatternP2P object.
202 */
203 void
205
206 /**
207 * @brief Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
208 */
209 void
211
212 /**
213 * @brief Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
214 */
215 void
217
218
219 /**
220 * @brief Computes the cell-level stiffness matrix.
221 */
222 void
223 computeCellStiffnessMatrix(const unsigned int quadratureID,
224 const unsigned int cellsBlockSize,
225 const bool basisType = false,
226 const bool ceoffType = true);
227
228 void
229 computeCellMassMatrix(const unsigned int quadratureID,
230 const unsigned int cellsBlockSize,
231 const bool basisType = false,
232 const bool ceoffType = true);
233
234 void
236 const std::pair<unsigned int, unsigned int> cellRangeTotal,
239 &weightedCellMassMatrix) const;
240
241 void
243 const std::pair<unsigned int, unsigned int> cellRangeTotal,
246 &weightedCellNjGradNiMatrix) const;
247
248 void
250 const std::pair<unsigned int, unsigned int> cellRangeTotal,
253 &weightedCellNjGradNiPlusNiGradNjMatrix) const;
254
255 void
256 computeInverseSqrtMassVector(const bool basisType = true,
257 const bool ceoffType = false);
258
259 /**
260 * @brief Computes the stiffness matrix
261 * \grad Ni . \grad Ni.
262 */
263 void
264 computeStiffnessVector(const bool basisType = true,
265 const bool ceoffType = false);
266
267 /**
268 * @brief Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in reinit.
269 */
270 void
271 resizeTempStorage(const bool isResizeTempStorageForInerpolation,
272 const bool isResizeTempStorageForCellMatrices);
273
274 /**
275 * @brief Number of quadrature points per cell for the quadratureID set in reinit.
276 */
277 unsigned int
279
280 /**
281 * @brief Number of vectors set in reinit.
282 */
283 unsigned int
284 nVectors() const;
285 /**
286 * @brief Number of DoFs per cell for the dofHandlerID set in init.
287 */
288 unsigned int
290
291 /**
292 * @brief Number of locally owned cells on the current processor.
293 */
294 unsigned int
295 nCells() const;
296
297 /**
298 * @brief Number of DoFs on the current processor, locally owned + ghosts.
299 */
300 unsigned int
302
303 /**
304 * @brief Number of locally owned DoFs on the current processor.
305 */
306 unsigned int
307 nOwnedDofs() const;
308
309 /**
310 * @brief Shape function values at quadrature points.
311 * @param[in] transpose if false the the data is indexed as [iQuad *
312 * d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
313 * d_nQuadsPerCell + iQuad].
314 */
316 shapeFunctionData(bool transpose = false) const;
317
318 /**
319 * @brief Shape function gradient values at quadrature points.
320 * @param[in] transpose if false the the data is indexed as [iDim *
321 * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
322 * if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
323 * iNode * d_nQuadsPerCell + iQuad].
324 */
326 shapeFunctionGradientData(bool transpose = false) const;
327
328 /**
329 * @brief Inverse Jacobian matrices, for cartesian cells returns the
330 * diagonal elements of the inverse Jacobian matrices for each cell, for
331 * affine cells returns the 3x3 inverse Jacobians for each cell otherwise
332 * returns the 3x3 inverse Jacobians at each quad point for each cell.
333 */
336
337 /**
338 * @brief determinant of Jacobian times the quadrature weight at each
339 * quad point for each cell.
340 */
342 JxW() const;
343
344 /**
345 * @brief quad point coordinates for each cell.
346 */
347 const dftfe::utils::MemoryStorage<ValueTypeBasisData,
349 quadPoints() const;
350
351 /**
352 * @brief Shape function values at quadrature points in ValueTypeBasisData.
353 * @param[in] transpose if false the the data is indexed as [iQuad *
354 * d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
355 * d_nQuadsPerCell + iQuad].
356 */
357 const auto &
358 shapeFunctionBasisData(bool transpose = false) const
359 {
360 if constexpr (std::is_same<ValueTypeBasisCoeff,
361 ValueTypeBasisData>::value)
362 {
363 return transpose ?
366 }
367 else
368 {
369 return transpose ?
371 ->second :
373 }
374 }
375 /**
376 * @brief Shape function gradient values at quadrature points in ValueTypeBasisData.
377 * @param[in] transpose if false the the data is indexed as [iDim *
378 * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
379 * if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
380 * iNode * d_nQuadsPerCell + iQuad].
381 */
382 const auto &
383 shapeFunctionGradientBasisData(bool transpose = false) const
384 {
385 if constexpr (std::is_same<ValueTypeBasisCoeff,
386 ValueTypeBasisData>::value)
387 {
388 return transpose ?
390 ->second :
392 }
393 else
394 {
395 return transpose ?
397 .find(d_quadratureID)
398 ->second :
400 ->second;
401 }
402 }
403
404 /**
405 * @brief Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the
406 * diagonal elements of the inverse Jacobian matrices for each cell, for
407 * affine cells returns the 3x3 inverse Jacobians for each cell otherwise
408 * returns the 3x3 inverse Jacobians at each quad point for each cell.
409 */
410 const auto &
412 {
413 if constexpr (std::is_same<ValueTypeBasisCoeff,
414 ValueTypeBasisData>::value)
415 {
418 ->second;
419 }
420 else
421 {
424 ->second;
425 }
426 }
427
428 /**
429 * @brief determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each
430 * quad point for each cell.
431 */
432 const auto &
434 {
435 if constexpr (std::is_same<ValueTypeBasisCoeff,
436 ValueTypeBasisData>::value)
437 {
438 return d_JxWData.find(d_quadratureID)->second;
439 }
440 else
441 {
442 return d_JxWBasisData.find(d_quadratureID)->second;
443 }
444 }
445
446 /**
447 * @brief Cell level stiffness matrix in ValueTypeBasisCoeff
448 */
449 const auto &
451 {
452 if constexpr (std::is_same<ValueTypeBasisCoeff,
453 ValueTypeBasisData>::value)
454 {
456 }
457 else
458 {
460 }
461 }
462
463
464 /**
465 * @brief Cell level stiffness matrix in ValueTypeBasisData
466 */
469
470
471 /**
472 * @brief Cell level mass matrix in ValueTypeBasisCoeff
473 */
474 const auto &
476 {
477 if constexpr (std::is_same<ValueTypeBasisCoeff,
478 ValueTypeBasisData>::value)
479 {
481 }
482 else
483 {
485 }
486 }
487
488
489 /**
490 * @brief Cell level mass matrix in ValueTypeBasisData
491 */
494
495
496 /**
497 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
498 */
499 const auto &
501 {
502 if constexpr (std::is_same<ValueTypeBasisCoeff,
503 ValueTypeBasisData>::value)
504 {
506 }
507 else
508 {
510 }
511 }
512
513 /**
514 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff
515 */
516 const auto &
518 {
519 if constexpr (std::is_same<ValueTypeBasisCoeff,
520 ValueTypeBasisData>::value)
521 {
523 }
524 else
525 {
527 }
528 }
529
530 /**
531 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff
532 */
533 const auto &
535 {
536 if constexpr (std::is_same<ValueTypeBasisCoeff,
537 ValueTypeBasisData>::value)
538 {
540 }
541 else
542 {
544 }
545 }
546
547 /**
548 * @brief Cell level diagonal mass matrix in ValueTypeBasisCoeff
549 */
550 const auto &
552 {
553 if constexpr (std::is_same<ValueTypeBasisCoeff,
554 ValueTypeBasisData>::value)
555 {
557 }
558 else
559 {
561 }
562 }
563
564
565 /**
566 * @brief Inverse sqrt diagonal 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 /**
585 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
586 */
587 const auto &
589 {
590 if constexpr (std::is_same<ValueTypeBasisCoeff,
591 ValueTypeBasisData>::value)
592 {
594 }
595 else
596 {
598 }
599 }
600
601
602 /**
603 * @brief sqrt diagonal mass matrix in ValueTypeBasisCoeff
604 */
605 const auto &
607 {
608 if constexpr (std::is_same<ValueTypeBasisCoeff,
609 ValueTypeBasisData>::value)
610 {
612 }
613 else
614 {
616 }
617 }
618
619 /**
620 * @brief diagonal mass matrix in ValueTypeBasisCoeff
621 */
622 const auto &
624 {
625 if constexpr (std::is_same<ValueTypeBasisCoeff,
626 ValueTypeBasisData>::value)
627 {
629 }
630 else
631 {
633 }
634 }
635
636
637 /**
638 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData
639 */
642
643
644 /**
645 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
646 */
649
650 /**
651 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
652 */
655 memorySpace> &
657
658 /**
659 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisData
660 */
663
664
665 /**
666 * @brief Cell level diagonal mass matrix in ValueTypeBasisData
667 */
670
671 /**
672 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisData
673 */
676
677 /**
678 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
679 */
682
683 /**
684 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
685 */
688 memorySpace> &
690
691 /**
692 * @brief sqrt diagonal mass matrix in ValueTypeBasisData
693 */
696
697 /**
698 * @brief diagonal mass matrix in ValueTypeBasisData
699 */
702
703 /**
704 * @brief diagonal stiffness matrix in ValueTypeBasisData
705 */
708
709 /**
710 * @brief diagonal inverse stiffness matrix in ValueTypeBasisData
711 */
714
715 /**
716 * @brief diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
717 */
720
721 /**
722 * @brief diagonal sqrt stiffness matrix in ValueTypeBasisData
723 */
726
727 /**
728 * @brief returns 2 if all cells on current processor are Cartesian,
729 * 1 if all cells on current processor are affine and 0 otherwise.
730 */
731 unsigned int
733
734 /**
735 * @brief returns the deal.ii cellID corresponing to given cell Index.
736 * @param[in] iElem cell Index
737 */
738 dealii::CellId
739 cellID(const unsigned int iElem) const;
740 /**
741 * @brief returns the deal.ii cell_iterator corresponing to given cell Index.
742 * @param[in] iElem cell Index
743 */
744
745 dealii::DoFHandler<3>::active_cell_iterator
746 getCellIterator(const unsigned int iElem) const;
747
748 /**
749 * @brief returns the cell index corresponding to given deal.ii cellID.
750 * @param[in] iElem cell Index
751 */
752 unsigned int
753 cellIndex(const dealii::CellId cellid) const;
754
755 /**
756 * @brief Creates a multivector.
757 * @param[in] blocksize Number of vectors in the multivector.
758 * @param[out] multiVector the created multivector.
759 */
760 void
762 const unsigned int blocksize,
764 &multiVector) const;
765
766 /**
767 * @brief Creates a multivector.
768 * @param[in] blocksize Number of vectors in the multivector.
769 * @param[out] multiVector the created multivector.
770 */
771 void
773 const unsigned int blocksize,
776 memorySpace> &multiVector) const;
777
778 /**
779 * @brief Creates scratch multivectors.
780 * @param[in] vecBlockSize Number of vectors in the multivector.
781 * @param[out] numMultiVecs number of scratch multivectors needed with
782 * this vecBlockSize.
783 */
784 void
785 createScratchMultiVectors(const unsigned int vecBlockSize,
786 const unsigned int numMultiVecs = 1) const;
787
788 /**
789 * @brief Creates single precision scratch multivectors.
790 * @param[in] vecBlockSize Number of vectors in the multivector.
791 * @param[out] numMultiVecs number of scratch multivectors needed with
792 * this vecBlockSize.
793 */
794 void
796 const unsigned int vecBlockSize,
797 const unsigned int numMultiVecs = 1) const;
798
799 /**
800 * @brief Clears scratch multivectors.
801 */
802 void
804
805 /**
806 * @brief Gets scratch multivectors.
807 * @param[in] vecBlockSize Number of vectors in the multivector.
808 * @param[out] numMultiVecs index of the multivector among those with the
809 * same vecBlockSize.
810 */
812 getMultiVector(const unsigned int vecBlockSize,
813 const unsigned int index = 0) const;
814
815 /**
816 * @brief Gets single precision scratch multivectors.
817 * @param[in] vecBlockSize Number of vectors in the multivector.
818 * @param[out] numMultiVecs index of the multivector among those with the
819 * same vecBlockSize.
820 */
823 memorySpace> &
824 getMultiVectorSinglePrec(const unsigned int vecBlockSize,
825 const unsigned int index = 0) const;
826
827 /**
828 * @brief Apply constraints on given multivector.
829 * @param[inout] multiVector the given multivector.
830 */
831 void
833 memorySpace> &multiVector,
834 unsigned int constraintIndex =
835 std::numeric_limits<unsigned int>::max()) const;
836
837
838
839 /**
840 * @brief Return the underlying deal.II matrixfree object.
841 */
842 const dealii::MatrixFree<3, ValueTypeBasisData> &
844
845 /**
846 * @brief Return the underlying deal.II dofhandler object.
847 */
848 const dealii::DoFHandler<3> &
850
851
852
853 std::vector<dftUtils::constraintMatrixInfo<memorySpace>> d_constraintInfo;
854 unsigned int d_nOMPThreads;
855 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
857 const dealii::MatrixFree<3, ValueTypeBasisData> *d_matrixFreeDataPtr;
861 std::map<unsigned int,
862 dftfe::utils::MemoryStorage<ValueTypeBasisData,
867 std::vector<dealii::CellId> d_cellIndexToCellIdMap;
868 std::vector<dealii::DoFHandler<3>::active_cell_iterator>
870 std::map<dealii::CellId, unsigned int> d_cellIdToCellIndexMap;
871 std::map<unsigned int,
874 std::map<unsigned int,
877 std::map<unsigned int,
880 std::map<unsigned int,
883 std::map<unsigned int,
886 std::map<unsigned int,
889 std::map<unsigned int,
892
893 std::map<unsigned int,
896 std::map<unsigned int,
899 std::map<unsigned int,
902 std::map<unsigned int,
905 std::map<unsigned int,
908 std::map<unsigned int,
911
922
925 memorySpace>
949 memorySpace>
961 mutable std::map<
962 unsigned int,
963 std::vector<
966
983
1000
1001
1002 mutable std::map<
1003 unsigned int,
1006 memorySpace>>>
1008
1009 std::vector<unsigned int> d_quadratureIDsVector;
1010 unsigned int d_quadratureID;
1011 unsigned int d_quadratureIndex;
1012 std::vector<unsigned int> d_nQuadsPerCell;
1013 unsigned int d_dofHandlerID;
1014 unsigned int d_nVectors;
1015 unsigned int d_nCells;
1016 unsigned int d_cellsBlockSize;
1017 unsigned int d_nDofsPerCell;
1018 unsigned int d_localSize;
1022 std::vector<UpdateFlags> d_updateFlags;
1023
1024 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
1026
1027
1028 /**
1029 * @brief Interpolate process level nodal data to cell level quadrature data.
1030 * @param[in] nodalData process level nodal data, the multivector should
1031 * already have ghost data and constraints should have been applied.
1032 * @param[out] quadratureValues Cell level quadrature values, indexed by
1033 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1034 * @param[out] quadratureGradients Cell level quadrature gradients,
1035 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1036 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1037 */
1038 void
1040 memorySpace> &nodalData,
1041 ValueTypeBasisCoeff *quadratureValues,
1042 ValueTypeBasisCoeff *quadratureGradients = NULL) const;
1043
1044 // FIXME Untested function
1045 /**
1046 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1047 * @param[in] quadratureValues Cell level quadrature values, indexed by
1048 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1049 * @param[in] quadratureGradients Cell level quadrature gradients,
1050 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1051 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1052 * @param[out] nodalData process level nodal data.
1053 */
1054 void
1056 ValueTypeBasisCoeff *quadratureValues,
1057 ValueTypeBasisCoeff *quadratureGradients,
1059 &nodalData,
1061 &mapQuadIdToProcId) const;
1062
1063 /**
1064 * @brief Get cell level nodal data from process level nodal data.
1065 * @param[in] nodalData process level nodal data, the multivector should
1066 * already have ghost data and constraints should have been applied.
1067 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1068 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1069 */
1070 void
1073 & nodalData,
1074 ValueTypeBasisCoeff *cellNodalDataPtr) const;
1075 // FIXME Untested function
1076 /**
1077 * @brief Accumulate cell level nodal data into process level nodal data.
1078 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1079 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1080 * @param[out] nodalData process level nodal data.
1081 */
1082 void
1084 const ValueTypeBasisCoeff *cellNodalDataPtr,
1086 &nodalData) const;
1087
1088 /**
1089 * @brief Interpolate process level nodal data to cell level quadrature data.
1090 * @param[in] nodalData process level nodal data, the multivector should
1091 * already have ghost data and constraints should have been applied.
1092 * @param[out] quadratureValues Cell level quadrature values, indexed by
1093 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1094 * @param[out] quadratureGradients Cell level quadrature gradients,
1095 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1096 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1097 * @param[in] cellRange the range of cells for which interpolation has to
1098 * be done.
1099 */
1100 void
1102 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1103 memorySpace> &nodalData,
1104 ValueTypeBasisCoeff * quadratureValues,
1105 ValueTypeBasisCoeff * quadratureGradients,
1106 const std::pair<unsigned int, unsigned int> cellRange) const;
1107
1108 /**
1109 * @brief Interpolate cell level nodal data to cell level quadrature data.
1110 * @param[in] nodalData cell level nodal data, the multivector should
1111 * already have ghost data and constraints should have been applied.
1112 * @param[out] quadratureValues Cell level quadrature values, indexed by
1113 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1114 * @param[out] quadratureGradients Cell level quadrature gradients,
1115 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1116 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1117 * @param[in] cellRange the range of cells for which interpolation has to
1118 * be done.
1119 */
1120 void
1122 const ValueTypeBasisCoeff * nodalData,
1123 ValueTypeBasisCoeff * quadratureValues,
1124 ValueTypeBasisCoeff * quadratureGradients,
1125 const std::pair<unsigned int, unsigned int> cellRange) const;
1126
1127 // FIXME Untested function
1128 /**
1129 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1130 * @param[in] quadratureValues Cell level quadrature values, indexed by
1131 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1132 * @param[in] quadratureGradients Cell level quadrature gradients,
1133 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1134 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1135 * @param[out] nodalData process level nodal data.
1136 * @param[in] cellRange the range of cells for which integration has to be
1137 * done.
1138 */
1139 void
1141 const ValueTypeBasisCoeff *quadratureValues,
1142 const ValueTypeBasisCoeff *quadratureGradients,
1144 &nodalData,
1146 & mapQuadIdToProcId,
1147 const std::pair<unsigned int, unsigned int> cellRange) const;
1148
1149
1150 /**
1151 * @brief Get cell level nodal data from process level nodal data.
1152 * @param[in] nodalData process level nodal data, the multivector should
1153 * already have ghost data and constraints should have been applied.
1154 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1155 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1156 * @param[in] cellRange the range of cells for which extraction has to be
1157 * done.
1158 */
1159 void
1161 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1162 memorySpace> &nodalData,
1163 ValueTypeBasisCoeff * cellNodalDataPtr,
1164 const std::pair<unsigned int, unsigned int> cellRange) const;
1165
1166 // FIXME Untested function
1167 /**
1168 * @brief Accumulate cell level nodal data into process level nodal data.
1169 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1170 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1171 * @param[out] nodalData process level nodal data.
1172 * @param[in] cellRange the range of cells for which extraction has to be
1173 * done.
1174 */
1175 void
1177 const ValueTypeBasisCoeff *cellNodalDataPtr,
1179 & nodalData,
1180 const std::pair<unsigned int, unsigned int> cellRange) const;
1181 };
1182 } // end of namespace basis
1183} // end of namespace dftfe
1184
1185#endif // dftfeBasisOperations_h
unsigned int d_quadratureIndex
Definition FEBasisOperations.h:1011
dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > d_flattenedCellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:866
const auto & massVector() const
diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:623
const dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > & quadPoints() const
quad point coordinates for each cell.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseSqrtMassVectorBasisData() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisData.
const auto & cellMassVector() const
Cell level diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:551
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseMassVectorCoeffType
Definition FEBasisOperations.h:952
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:972
std::vector< unsigned int > d_quadratureIDsVector
Definition FEBasisOperations.h:1009
const auto & inverseSqrtMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:569
unsigned int nQuadsPerCell() const
Number of quadrature points per cell for the quadratureID set in reinit.
void computeCellMassMatrix(const unsigned int quadratureID, const unsigned int cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
bool areAllCellsAffine
Definition FEBasisOperations.h:1020
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:956
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseStiffnessVectorBasisData() const
diagonal inverse stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< dftfe::global_size_type, dftfe::utils::MemorySpace::HOST > d_cellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:860
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:993
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_cellInverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:926
unsigned int cellsTypeFlag() const
returns 2 if all cells on current processor are Cartesian, 1 if all cells on current processor are af...
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtMassVectorCoeffType
Definition FEBasisOperations.h:940
const dealii::MatrixFree< 3, ValueTypeBasisData > * d_matrixFreeDataPtr
Definition FEBasisOperations.h:857
const auto & cellInverseMassVector() const
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:517
unsigned int d_nDofsPerCell
Definition FEBasisOperations.h:1017
void initializeMPIPattern()
Constructs the MPIPatternP2P object.
unsigned int d_nVectors
Definition FEBasisOperations.h:1014
void interpolate(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients=NULL) const
Interpolate process level nodal data to cell level quadrature data.
std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > d_cellIndexToCellIteratorMap
Definition FEBasisOperations.h:869
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock
Definition FEBasisOperations.h:90
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseMassVectorBasisType
Definition FEBasisOperations.h:921
void reinit(const unsigned int &vecBlockSize, const unsigned int &cellBlockSize, const unsigned int &quadratureID, const bool isResizeTempStorageForInerpolation=true, const bool isResizeTempStorageForCellMatrices=false)
sets internal variables and optionally resizes internal temp storage for interpolation operations
const dealii::DoFHandler< 3 > & getDofHandler() const
Return the underlying deal.II dofhandler object.
std::map< dealii::CellId, unsigned int > d_cellIdToCellIndexMap
Definition FEBasisOperations.h:870
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:987
dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > & getMultiVectorSinglePrec(const unsigned int vecBlockSize, const unsigned int index=0) const
Gets single precision scratch multivectors.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_massVectorCoeffType
Definition FEBasisOperations.h:944
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellStiffnessMatrixBasisData() const
Cell level stiffness matrix in ValueTypeBasisData.
void computeWeightedCellMassMatrix(const std::pair< unsigned int, unsigned int > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellMassMatrix) const
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellInverseMassVectorBasisData() const
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
void initializeConstraints()
Initializes the constraintMatrixInfo object.
void integrateWithBasisKernel(const ValueTypeBasisCoeff *quadratureValues, const ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > &mapQuadIdToProcId, const std::pair< unsigned int, unsigned int > cellRange) const
Integrate cell level quadrature data times shape functions to process level nodal data.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_JxWData
Definition FEBasisOperations.h:876
~FEBasisOperations()=default
Default Destructor.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseStiffnessVectorBasisType
Definition FEBasisOperations.h:970
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:976
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientData
Definition FEBasisOperations.h:885
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:974
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & inverseJacobians() const
Inverse Jacobian matrices, for cartesian cells returns the diagonal elements of the inverse Jacobian ...
void extractToCellNodalData(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr) const
Get cell level nodal data from process level nodal data.
unsigned int d_quadratureID
Definition FEBasisOperations.h:1010
const auto & inverseJacobiansBasisData() const
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of...
Definition FEBasisOperations.h:411
const auto & cellStiffnessMatrix() const
Cell level stiffness matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:450
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
Definition FEBasisOperations.h:97
void resizeTempStorage(const bool isResizeTempStorageForInerpolation, const bool isResizeTempStorageForCellMatrices)
Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in re...
dftfe::utils::MemoryStorage< dftfe::global_size_type, dftfe::utils::MemorySpace::HOST > & getFlattenedMapsHost()
const auto & shapeFunctionBasisData(bool transpose=false) const
Shape function values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:358
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_stiffnessVectorBasisType
Definition FEBasisOperations.h:982
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock2
Definition FEBasisOperations.h:90
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:954
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > mpiPatternP2P
Definition FEBasisOperations.h:1025
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:995
void computeInverseSqrtMassVector(const bool basisType=true, const bool ceoffType=false)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:985
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:991
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseStiffnessVectorBasisType
Definition FEBasisOperations.h:980
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessVectorBasisType
Definition FEBasisOperations.h:968
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtMassVectorCoeffType
Definition FEBasisOperations.h:960
void accumulateFromCellNodalData(const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData) const
Accumulate cell level nodal data into process level nodal data.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & massVectorBasisData() const
diagonal mass matrix in ValueTypeBasisData
void reinitializeConstraints(std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector)
Reinitializes the constraintMatrixInfo object.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassMatrixCoeffType
Definition FEBasisOperations.h:919
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtMassVectorBasisType
Definition FEBasisOperations.h:958
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & shapeFunctionGradientData(bool transpose=false) const
Shape function gradient values at quadrature points.
std::vector< dealii::CellId > d_cellIndexToCellIdMap
Definition FEBasisOperations.h:867
void computeWeightedCellNjGradNiMatrix(const std::pair< unsigned int, unsigned int > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiMatrix) const
void clearScratchMultiVectors() const
Clears scratch multivectors.
std::vector< UpdateFlags > d_updateFlags
Definition FEBasisOperations.h:1022
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_inverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:950
void interpolateKernel(const ValueTypeBasisCoeff *nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< unsigned int, unsigned int > cellRange) const
Interpolate cell level nodal data to cell level quadrature data.
unsigned int nVectors() const
Number of vectors set in reinit.
const auto & sqrtMassVector() const
sqrt diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:606
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:999
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & JxW() const
determinant of Jacobian times the quadrature weight at each quad point for each cell.
unsigned int nCells() const
Number of locally owned cells on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessMatrixBasisType
Definition FEBasisOperations.h:913
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & cellInverseMassVectorBasisDataSinglePrec() const
Cell level inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessVectorCoeffType
Definition FEBasisOperations.h:989
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:930
unsigned int cellIndex(const dealii::CellId cellid) const
returns the cell index corresponding to given deal.ii cellID.
void clear()
Clears the FEBasisOperations internal storage.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisDataTranspose
Definition FEBasisOperations.h:910
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > > d_quadPoints
Definition FEBasisOperations.h:864
void distribute(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector, unsigned int constraintIndex=std::numeric_limits< unsigned int >::max()) const
Apply constraints on given multivector.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellSqrtMassVectorBasisData() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseMassVectorBasisType
Definition FEBasisOperations.h:946
const auto & cellSqrtMassVector() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:534
unsigned int d_dofHandlerID
Definition FEBasisOperations.h:1013
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisData
Definition FEBasisOperations.h:904
void computeStiffnessVector(const bool basisType=true, const bool ceoffType=false)
Computes the stiffness matrix \grad Ni . \grad Ni.
unsigned int d_locallyOwnedSize
Definition FEBasisOperations.h:1019
void computeWeightedCellNjGradNiPlusNiGradNjMatrix(const std::pair< unsigned int, unsigned int > cellRangeTotal, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weights, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > &weightedCellNjGradNiPlusNiGradNjMatrix) const
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsData
Definition FEBasisOperations.h:87
dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > zeroIndexVec
Definition FEBasisOperations.h:93
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataInternalLayout
Definition FEBasisOperations.h:882
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellValuesBlockCoeff
Definition FEBasisOperations.h:95
void accumulateFromCellNodalDataKernel(const ValueTypeBasisCoeff *cellNodalDataPtr, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, const std::pair< unsigned int, unsigned int > cellRange) const
Accumulate cell level nodal data into process level nodal data.
unsigned int nRelaventDofs() const
Number of DoFs on the current processor, locally owned + ghosts.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellStiffnessMatrixCoeffType
Definition FEBasisOperations.h:915
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionData
Definition FEBasisOperations.h:879
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_stiffnessVectorCoeffType
Definition FEBasisOperations.h:997
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellMassVectorBasisData() const
Cell level diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:932
std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > * d_constraintsVector
Definition FEBasisOperations.h:856
const auto & inverseMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:588
void init(dealii::MatrixFree< 3, ValueTypeBasisData > &matrixFreeData, std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > &constraintsVector, const unsigned int &dofHandlerID, const std::vector< unsigned int > &quadratureID, const std::vector< UpdateFlags > updateFlags)
fills required data structures for the given dofHandlerID
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_JxWBasisData
Definition FEBasisOperations.h:898
unsigned int d_cellsBlockSize
Definition FEBasisOperations.h:1016
void computeCellStiffnessMatrix(const unsigned int quadratureID, const unsigned int cellsBlockSize, const bool basisType=false, const bool ceoffType=true)
Computes the cell-level stiffness matrix.
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > & inverseMassVectorBasisDataSinglePrec() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellValuesBlock
Definition FEBasisOperations.h:90
dealii::DoFHandler< 3 >::active_cell_iterator getCellIterator(const unsigned int iElem) const
returns the deal.ii cell_iterator corresponing to given cell Index.
void initializeShapeFunctionAndJacobianData()
Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtMassVectorBasisType
Definition FEBasisOperations.h:938
void initializeShapeFunctionAndJacobianBasisData()
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseMassVectorCoeffType
Definition FEBasisOperations.h:928
unsigned int nDofsPerCell() const
Number of DoFs per cell for the dofHandlerID set in init.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_sqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:978
const auto & JxWBasisData() const
determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each quad point for each...
Definition FEBasisOperations.h:433
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisData
Definition FEBasisOperations.h:901
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtStiffnessVectorBasisData() const
diagonal sqrt stiffness matrix in ValueTypeBasisData
void init(const FEBasisOperations< ValueTypeBasisCoeff, ValueTypeBasisData, memorySpaceSrc > &basisOperationsSrc)
fills required data structures from another FEBasisOperations object
std::map< unsigned int, std::vector< dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > > > scratchMultiVectorsSinglePrec
Definition FEBasisOperations.h:1007
bool areAllCellsCartesian
Definition FEBasisOperations.h:1021
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellNodalData
Definition FEBasisOperations.h:87
void createScratchMultiVectors(const unsigned int vecBlockSize, const unsigned int numMultiVecs=1) const
Creates scratch multivectors.
unsigned int d_nCells
Definition FEBasisOperations.h:1015
unsigned int d_localSize
Definition FEBasisOperations.h:1018
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassVectorBasisType
Definition FEBasisOperations.h:934
void initializeFlattenedIndexMaps()
Initializes indexset maps from process level indices to cell level indices for multivectors.
const auto & cellInverseSqrtMassVector() const
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:500
void initializeIndexMaps()
Initializes indexset maps from process level indices to cell level indices for a single vector,...
std::vector< dftUtils::constraintMatrixInfo< memorySpace > > d_constraintInfo
Definition FEBasisOperations.h:853
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseSqrtStiffnessVectorBasisData() const
diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionDataTranspose
Definition FEBasisOperations.h:888
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_inverseJacobianData
Definition FEBasisOperations.h:873
FEBasisOperations(std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr)
Constructor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellMassVectorCoeffType
Definition FEBasisOperations.h:936
dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > & getMultiVector(const unsigned int vecBlockSize, const unsigned int index=0) const
Gets scratch multivectors.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseMassVectorBasisData() const
Inverse diagonal mass matrix in ValueTypeBasisData.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellMatrixBlock
Definition FEBasisOperations.h:91
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_inverseJacobianBasisData
Definition FEBasisOperations.h:895
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & sqrtMassVectorBasisData() const
sqrt diagonal mass matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassMatrixBasisType
Definition FEBasisOperations.h:917
void createMultiVectorSinglePrec(const unsigned int blocksize, dftfe::linearAlgebra::MultiVector< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisCoeff >::type, memorySpace > &multiVector) const
Creates a multivector.
void createScratchMultiVectorsSinglePrec(const unsigned int vecBlockSize, const unsigned int numMultiVecs=1) const
Creates single precision scratch multivectors.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & stiffnessVectorBasisData() const
diagonal stiffness matrix in ValueTypeBasisData
dealii::CellId cellID(const unsigned int iElem) const
returns the deal.ii cellID corresponing to given cell Index.
std::map< unsigned int, std::vector< dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > > > scratchMultiVectors
Definition FEBasisOperations.h:965
void createMultiVector(const unsigned int blocksize, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector) const
Creates a multivector.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsDataNonAffine
Definition FEBasisOperations.h:88
const dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > & shapeFunctionData(bool transpose=false) const
Shape function values at quadrature points.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataTranspose
Definition FEBasisOperations.h:891
std::vector< unsigned int > d_nQuadsPerCell
Definition FEBasisOperations.h:1012
const auto & cellMassMatrix() const
Cell level mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:475
const auto & shapeFunctionGradientBasisData(bool transpose=false) const
Shape function gradient values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:383
unsigned int d_nOMPThreads
Definition FEBasisOperations.h:854
const dealii::MatrixFree< 3, ValueTypeBasisData > & matrixFreeData() const
Return the underlying deal.II matrixfree object.
void extractToCellNodalDataKernel(const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *cellNodalDataPtr, const std::pair< unsigned int, unsigned int > cellRange) const
Get cell level nodal data from process level nodal data.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellInverseSqrtMassVectorBasisData() const
Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & cellMassMatrixBasisData() const
Cell level mass matrix in ValueTypeBasisData.
void integrateWithBasis(ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > &mapQuadIdToProcId) const
Integrate cell level quadrature data times shape functions to process level nodal data.
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_massVectorBasisType
Definition FEBasisOperations.h:942
void interpolateKernel(const dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &nodalData, ValueTypeBasisCoeff *quadratureValues, ValueTypeBasisCoeff *quadratureGradients, const std::pair< unsigned int, unsigned int > cellRange) const
Interpolate process level nodal data to cell level quadrature data.
std::map< unsigned int, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisDataTranspose
Definition FEBasisOperations.h:907
unsigned int nOwnedDofs() const
Number of locally owned DoFs on the current processor.
Definition BLASWrapper.h:35
An class template to encapsulate a MultiVector. A MultiVector is a collection of vectors belonging t...
Definition MultiVector.h:127
Definition MemoryStorage.h:33
Definition FEBasisOperations.h:30
UpdateFlags operator|(const UpdateFlags f1, const UpdateFlags f2)
Definition FEBasisOperations.h:49
UpdateFlags & operator|=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:58
UpdateFlags
Definition FEBasisOperations.h:32
@ update_values
Definition FEBasisOperations.h:35
@ update_quadpoints
Definition FEBasisOperations.h:41
@ update_transpose
Definition FEBasisOperations.h:39
@ update_inversejacobians
Definition FEBasisOperations.h:43
@ update_jxw
Definition FEBasisOperations.h:45
@ update_default
Definition FEBasisOperations.h:33
@ update_gradients
Definition FEBasisOperations.h:37
UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
Definition FEBasisOperations.h:65
UpdateFlags & operator&=(UpdateFlags &f1, const UpdateFlags f2)
Definition FEBasisOperations.h:73
MemorySpace
Definition MemorySpaceType.h:33
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
unsigned long int global_size_type
Definition TypeConfig.h:7
T type
Definition dftfeDataTypes.h:115