DFT-FE 1.3.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
FEBasisOperations.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17
18#ifndef dftfeFEBasisOperations_h
19#define dftfeFEBasisOperations_h
20
21#include <MultiVector.h>
22#include <headers.h>
24#include <DeviceTypeConfig.h>
25#include <BLASWrapper.h>
26
27namespace dftfe
28{
29 namespace basis
30 {
32 {
34
35 update_values = 0x0001,
36
38
40
42
44
45 update_jxw = 0x0020,
46 };
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 quad point coordinates for each cell.
368 */
369 const dftfe::utils::MemoryStorage<ValueTypeBasisData,
372
373 /**
374 * @brief Shape function values at quadrature points in ValueTypeBasisData.
375 * @param[in] transpose if false the the data is indexed as [iQuad *
376 * d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
377 * d_nQuadsPerCell + iQuad].
378 */
379 const auto &
380 shapeFunctionBasisData(bool transpose = false) const
381 {
382 if constexpr (std::is_same<ValueTypeBasisCoeff,
383 ValueTypeBasisData>::value)
384 {
385 return transpose ?
388 }
389 else
390 {
391 return transpose ?
393 ->second :
395 }
396 }
397 /**
398 * @brief Shape function gradient values at quadrature points in ValueTypeBasisData.
399 * @param[in] transpose if false the the data is indexed as [iDim *
400 * d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
401 * if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
402 * iNode * d_nQuadsPerCell + iQuad].
403 */
404 const auto &
405 shapeFunctionGradientBasisData(bool transpose = false) const
406 {
407 if constexpr (std::is_same<ValueTypeBasisCoeff,
408 ValueTypeBasisData>::value)
409 {
410 return transpose ?
412 ->second :
414 }
415 else
416 {
417 return transpose ?
419 .find(d_quadratureID)
420 ->second :
422 ->second;
423 }
424 }
425
426 /**
427 * @brief Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the
428 * diagonal elements of the inverse Jacobian matrices for each cell, for
429 * affine cells returns the 3x3 inverse Jacobians for each cell otherwise
430 * returns the 3x3 inverse Jacobians at each quad point for each cell.
431 */
432 const auto &
434 {
435 if constexpr (std::is_same<ValueTypeBasisCoeff,
436 ValueTypeBasisData>::value)
437 {
440 ->second;
441 }
442 else
443 {
446 ->second;
447 }
448 }
449
450 /**
451 * @brief determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each
452 * quad point for each cell.
453 */
454 const auto &
456 {
457 if constexpr (std::is_same<ValueTypeBasisCoeff,
458 ValueTypeBasisData>::value)
459 {
460 return d_JxWData.find(d_quadratureID)->second;
461 }
462 else
463 {
464 return d_JxWBasisData.find(d_quadratureID)->second;
465 }
466 }
467
468 /**
469 * @brief Cell level stiffness matrix in ValueTypeBasisCoeff
470 */
471 const auto &
473 {
474 if constexpr (std::is_same<ValueTypeBasisCoeff,
475 ValueTypeBasisData>::value)
476 {
478 }
479 else
480 {
482 }
483 }
484
485
486 /**
487 * @brief Cell level stiffness matrix in ValueTypeBasisData
488 */
491
492
493 /**
494 * @brief Cell level mass matrix in ValueTypeBasisCoeff
495 */
496 const auto &
498 {
499 if constexpr (std::is_same<ValueTypeBasisCoeff,
500 ValueTypeBasisData>::value)
501 {
503 }
504 else
505 {
507 }
508 }
509
510
511 /**
512 * @brief Cell level mass matrix in ValueTypeBasisData
513 */
516
517
518 /**
519 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
520 */
521 const auto &
523 {
524 if constexpr (std::is_same<ValueTypeBasisCoeff,
525 ValueTypeBasisData>::value)
526 {
528 }
529 else
530 {
532 }
533 }
534
535 /**
536 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff
537 */
538 const auto &
540 {
541 if constexpr (std::is_same<ValueTypeBasisCoeff,
542 ValueTypeBasisData>::value)
543 {
545 }
546 else
547 {
549 }
550 }
551
552 /**
553 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff
554 */
555 const auto &
557 {
558 if constexpr (std::is_same<ValueTypeBasisCoeff,
559 ValueTypeBasisData>::value)
560 {
562 }
563 else
564 {
566 }
567 }
568
569 /**
570 * @brief Cell level diagonal mass matrix in ValueTypeBasisCoeff
571 */
572 const auto &
574 {
575 if constexpr (std::is_same<ValueTypeBasisCoeff,
576 ValueTypeBasisData>::value)
577 {
579 }
580 else
581 {
583 }
584 }
585
586
587 /**
588 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
589 */
590 const auto &
592 {
593 if constexpr (std::is_same<ValueTypeBasisCoeff,
594 ValueTypeBasisData>::value)
595 {
597 }
598 else
599 {
601 }
602 }
603
604
605
606 /**
607 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
608 */
609 const auto &
611 {
612 if constexpr (std::is_same<ValueTypeBasisCoeff,
613 ValueTypeBasisData>::value)
614 {
616 }
617 else
618 {
620 }
621 }
622
623
624 /**
625 * @brief sqrt diagonal mass matrix in ValueTypeBasisCoeff
626 */
627 const auto &
629 {
630 if constexpr (std::is_same<ValueTypeBasisCoeff,
631 ValueTypeBasisData>::value)
632 {
634 }
635 else
636 {
638 }
639 }
640
641 /**
642 * @brief diagonal mass matrix in ValueTypeBasisCoeff
643 */
644 const auto &
646 {
647 if constexpr (std::is_same<ValueTypeBasisCoeff,
648 ValueTypeBasisData>::value)
649 {
651 }
652 else
653 {
655 }
656 }
657
658
659 /**
660 * @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData
661 */
664
665
666 /**
667 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
668 */
671
672 /**
673 * @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
674 */
677 memorySpace> &
679
680 /**
681 * @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisData
682 */
685
686
687 /**
688 * @brief Cell level diagonal mass matrix in ValueTypeBasisData
689 */
692
693 /**
694 * @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisData
695 */
698
699 /**
700 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
701 */
704
705 /**
706 * @brief Inverse diagonal mass matrix in ValueTypeBasisData
707 */
710 memorySpace> &
712
713 /**
714 * @brief sqrt diagonal mass matrix in ValueTypeBasisData
715 */
718
719 /**
720 * @brief diagonal mass matrix in ValueTypeBasisData
721 */
724
725 /**
726 * @brief diagonal stiffness matrix in ValueTypeBasisData
727 */
730
731 /**
732 * @brief diagonal inverse stiffness matrix in ValueTypeBasisData
733 */
736
737 /**
738 * @brief diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
739 */
742
743 /**
744 * @brief diagonal sqrt stiffness matrix in ValueTypeBasisData
745 */
748
749 /**
750 * @brief returns 2 if all cells on current processor are Cartesian,
751 * 1 if all cells on current processor are affine and 0 otherwise.
752 */
755
756 /**
757 * @brief returns the deal.ii cellID corresponing to given cell Index.
758 * @param[in] iElem cell Index
759 */
760 dealii::CellId
761 cellID(const dftfe::uInt iElem) const;
762 /**
763 * @brief returns the deal.ii cell_iterator corresponing to given cell Index.
764 * @param[in] iElem cell Index
765 */
766
767 dealii::DoFHandler<3>::active_cell_iterator
768 getCellIterator(const dftfe::uInt iElem) const;
769
770 /**
771 * @brief returns the cell index corresponding to given deal.ii cellID.
772 * @param[in] iElem cell Index
773 */
775 cellIndex(const dealii::CellId cellid) const;
776
777 /**
778 * @brief Creates a multivector.
779 * @param[in] blocksize Number of vectors in the multivector.
780 * @param[out] multiVector the created multivector.
781 */
782 void
784 const dftfe::uInt blocksize,
786 &multiVector) const;
787
788 /**
789 * @brief Creates a multivector.
790 * @param[in] blocksize Number of vectors in the multivector.
791 * @param[out] multiVector the created multivector.
792 */
793 void
795 const dftfe::uInt blocksize,
798 memorySpace> &multiVector) const;
799
800 /**
801 * @brief Creates scratch multivectors.
802 * @param[in] vecBlockSize Number of vectors in the multivector.
803 * @param[out] numMultiVecs number of scratch multivectors needed with
804 * this vecBlockSize.
805 */
806 void
808 const dftfe::uInt numMultiVecs = 1) const;
809
810 /**
811 * @brief Creates single precision scratch multivectors.
812 * @param[in] vecBlockSize Number of vectors in the multivector.
813 * @param[out] numMultiVecs number of scratch multivectors needed with
814 * this vecBlockSize.
815 */
816 void
818 const dftfe::uInt vecBlockSize,
819 const dftfe::uInt numMultiVecs = 1) const;
820
821 /**
822 * @brief Clears scratch multivectors.
823 */
824 void
826
827 /**
828 * @brief Gets scratch multivectors.
829 * @param[in] vecBlockSize Number of vectors in the multivector.
830 * @param[out] numMultiVecs index of the multivector among those with the
831 * same vecBlockSize.
832 */
834 getMultiVector(const dftfe::uInt vecBlockSize,
835 const dftfe::uInt index = 0) const;
836
837 /**
838 * @brief Gets single precision scratch multivectors.
839 * @param[in] vecBlockSize Number of vectors in the multivector.
840 * @param[out] numMultiVecs index of the multivector among those with the
841 * same vecBlockSize.
842 */
845 memorySpace> &
847 const dftfe::uInt index = 0) const;
848
849 /**
850 * @brief Apply constraints on given multivector.
851 * @param[inout] multiVector the given multivector.
852 */
853 void
855 memorySpace> &multiVector,
856 dftfe::uInt constraintIndex =
857 std::numeric_limits<dftfe::uInt>::max()) const;
858
859
860
861 /**
862 * @brief Return the underlying deal.II matrixfree object.
863 */
864 const dealii::MatrixFree<3, ValueTypeBasisData> &
866
867 /**
868 * @brief Return the underlying deal.II dofhandler object.
869 */
870 const dealii::DoFHandler<3> &
872
873
874
875 std::vector<dftUtils::constraintMatrixInfo<memorySpace>> d_constraintInfo;
877 std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
879 const dealii::MatrixFree<3, ValueTypeBasisData> *d_matrixFreeDataPtr;
882 std::map<dftfe::uInt,
883 dftfe::utils::MemoryStorage<ValueTypeBasisData,
886 dftfe::utils::MemoryStorage<ValueTypeBasisData,
891 std::vector<dealii::CellId> d_cellIndexToCellIdMap;
892 std::vector<dealii::DoFHandler<3>::active_cell_iterator>
894 std::map<dealii::CellId, dftfe::uInt> d_cellIdToCellIndexMap;
895 std::map<dftfe::uInt,
898 std::map<dftfe::uInt,
901 std::map<dftfe::uInt,
904 std::map<dftfe::uInt,
907 std::map<dftfe::uInt,
910 std::map<dftfe::uInt,
913 std::map<dftfe::uInt,
916
917 std::map<dftfe::uInt,
920 std::map<dftfe::uInt,
923 std::map<dftfe::uInt,
926 std::map<dftfe::uInt,
929 std::map<dftfe::uInt,
932 std::map<dftfe::uInt,
935
946
949 memorySpace>
973 memorySpace>
985 mutable std::map<
987 std::vector<
990
1007
1024
1025
1026 mutable std::map<
1030 memorySpace>>>
1032
1033 std::vector<dftfe::uInt> d_quadratureIDsVector;
1036 std::vector<dftfe::uInt> d_nQuadsPerCell;
1046 std::vector<UpdateFlags> d_updateFlags;
1047
1048 std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
1050
1051
1052 void
1054 distributedCPUVec<double> &nodalField,
1055 const dftfe::uInt dofHandlerId,
1056 const dftfe::uInt quadratureId,
1058 &quadratureValueData,
1060 &quadratureGradValueData,
1062 &quadratureHessianValueData,
1063 const bool isEvaluateGradData = false,
1064 const bool isEvaluateHessianData = false) const;
1065
1066 /**
1067 * @brief Interpolate process level nodal data to cell level quadrature data.
1068 * @param[in] nodalData process level nodal data, the multivector should
1069 * already have ghost data and constraints should have been applied.
1070 * @param[out] quadratureValues Cell level quadrature values, indexed by
1071 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1072 * @param[out] quadratureGradients Cell level quadrature gradients,
1073 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1074 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1075 */
1076 void
1078 memorySpace> &nodalData,
1079 ValueTypeBasisCoeff *quadratureValues,
1080 ValueTypeBasisCoeff *quadratureGradients = NULL) const;
1081
1082 // FIXME Untested function
1083 /**
1084 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1085 * @param[in] quadratureValues Cell level quadrature values, indexed by
1086 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1087 * @param[in] quadratureGradients Cell level quadrature gradients,
1088 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1089 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1090 * @param[out] nodalData process level nodal data.
1091 */
1092 void
1094 ValueTypeBasisCoeff *quadratureValues,
1095 ValueTypeBasisCoeff *quadratureGradients,
1097 &nodalData,
1099 &mapQuadIdToProcId) const;
1100
1101 /**
1102 * @brief Get cell level nodal data from process level nodal data.
1103 * @param[in] nodalData process level nodal data, the multivector should
1104 * already have ghost data and constraints should have been applied.
1105 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1106 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1107 */
1108 void
1111 &nodalData,
1112 ValueTypeBasisCoeff *cellNodalDataPtr) const;
1113 // FIXME Untested function
1114 /**
1115 * @brief Accumulate cell level nodal data into process level nodal data.
1116 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1117 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1118 * @param[out] nodalData process level nodal data.
1119 */
1120 void
1122 const ValueTypeBasisCoeff *cellNodalDataPtr,
1124 &nodalData) const;
1125
1126 /**
1127 * @brief Interpolate process level nodal data to cell level quadrature data.
1128 * @param[in] nodalData process level nodal data, the multivector should
1129 * already have ghost data and constraints should have been applied.
1130 * @param[out] quadratureValues Cell level quadrature values, indexed by
1131 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1132 * @param[out] 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[in] cellRange the range of cells for which interpolation has to
1136 * be done.
1137 */
1138 void
1140 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1141 memorySpace> &nodalData,
1142 ValueTypeBasisCoeff *quadratureValues,
1143 ValueTypeBasisCoeff *quadratureGradients,
1144 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1145
1146 /**
1147 * @brief Interpolate cell level nodal data to cell level quadrature data.
1148 * @param[in] nodalData cell level nodal data, the multivector should
1149 * already have ghost data and constraints should have been applied.
1150 * @param[out] quadratureValues Cell level quadrature values, indexed by
1151 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1152 * @param[out] quadratureGradients Cell level quadrature gradients,
1153 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1154 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1155 * @param[in] cellRange the range of cells for which interpolation has to
1156 * be done.
1157 */
1158 void
1160 const ValueTypeBasisCoeff *nodalData,
1161 ValueTypeBasisCoeff *quadratureValues,
1162 ValueTypeBasisCoeff *quadratureGradients,
1163 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1164
1165 // FIXME Untested function
1166 /**
1167 * @brief Integrate cell level quadrature data times shape functions to process level nodal data.
1168 * @param[in] quadratureValues Cell level quadrature values, indexed by
1169 * [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1170 * @param[in] quadratureGradients Cell level quadrature gradients,
1171 * indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
1172 * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
1173 * @param[out] nodalData process level nodal data.
1174 * @param[in] cellRange the range of cells for which integration has to be
1175 * done.
1176 */
1177 void
1179 const ValueTypeBasisCoeff *quadratureValues,
1180 const ValueTypeBasisCoeff *quadratureGradients,
1182 &nodalData,
1184 &mapQuadIdToProcId,
1185 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1186
1187
1188 /**
1189 * @brief Get cell level nodal data from process level nodal data.
1190 * @param[in] nodalData process level nodal data, the multivector should
1191 * already have ghost data and constraints should have been applied.
1192 * @param[out] cellNodalDataPtr Cell level nodal values, indexed by
1193 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1194 * @param[in] cellRange the range of cells for which extraction has to be
1195 * done.
1196 */
1197 void
1199 const dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
1200 memorySpace> &nodalData,
1201 ValueTypeBasisCoeff *cellNodalDataPtr,
1202 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1203
1204 // FIXME Untested function
1205 /**
1206 * @brief Accumulate cell level nodal data into process level nodal data.
1207 * @param[in] cellNodalDataPtr Cell level nodal values, indexed by
1208 * [iCell * d_nDofsPerCell * d_nVectors + iDoF * d_nVectors + iVec].
1209 * @param[out] nodalData process level nodal data.
1210 * @param[in] cellRange the range of cells for which extraction has to be
1211 * done.
1212 */
1213 void
1215 const ValueTypeBasisCoeff *cellNodalDataPtr,
1217 &nodalData,
1218 const std::pair<dftfe::uInt, dftfe::uInt> cellRange) const;
1219 };
1220 } // end of namespace basis
1221} // end of namespace dftfe
1222
1223#endif // dftfeBasisOperations_h
dftfe::uInt d_quadratureIndex
Definition FEBasisOperations.h:1035
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:1033
const auto & massVector() const
diagonal mass matrix in ValueTypeBasisCoeff
Definition FEBasisOperations.h:645
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:573
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseMassVectorCoeffType
Definition FEBasisOperations.h:976
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:996
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionBasisData
Definition FEBasisOperations.h:925
const auto & inverseSqrtMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:591
bool areAllCellsAffine
Definition FEBasisOperations.h:1044
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtMassVectorCoeffType
Definition FEBasisOperations.h:980
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & inverseStiffnessVectorBasisData() const
diagonal inverse stiffness matrix in ValueTypeBasisData
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_inverseSqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1017
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_cellInverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:950
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellSqrtMassVectorCoeffType
Definition FEBasisOperations.h:964
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_inverseJacobianData
Definition FEBasisOperations.h:897
const dealii::MatrixFree< 3, ValueTypeBasisData > * d_matrixFreeDataPtr
Definition FEBasisOperations.h:879
const auto & cellInverseMassVector() const
Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:539
void initializeMPIPattern()
Constructs the MPIPatternP2P object.
dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > d_cellCentroids
Definition FEBasisOperations.h:888
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:893
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > tempCellGradientsBlock
Definition FEBasisOperations.h:91
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseMassVectorBasisType
Definition FEBasisOperations.h:945
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:1043
std::map< dealii::CellId, dftfe::uInt > d_cellIdToCellIndexMap
Definition FEBasisOperations.h:894
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > > d_quadPoints
Definition FEBasisOperations.h:885
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:1011
dftfe::uInt nOwnedDofs() const
Number of locally owned DoFs on the current processor.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_massVectorCoeffType
Definition FEBasisOperations.h:968
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionGradientDataInternalLayout
Definition FEBasisOperations.h:906
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:915
~FEBasisOperations()=default
Default Destructor.
dftfe::uInt d_dofHandlerID
Definition FEBasisOperations.h:1037
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseStiffnessVectorBasisType
Definition FEBasisOperations.h:994
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_shapeFunctionData
Definition FEBasisOperations.h:903
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:1000
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtStiffnessVectorBasisType
Definition FEBasisOperations.h:998
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:912
const auto & inverseJacobiansBasisData() const
Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the diagonal elements of...
Definition FEBasisOperations.h:433
const auto & cellStiffnessMatrix() const
Cell level stiffness matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:472
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:1039
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:380
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:1006
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:978
std::shared_ptr< const utils::mpi::MPIPatternP2P< memorySpace > > mpiPatternP2P
Definition FEBasisOperations.h:1049
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtStiffnessVectorCoeffType
Definition FEBasisOperations.h:1019
void computeInverseSqrtMassVector(const bool basisType=true, const bool ceoffType=false)
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseStiffnessVectorCoeffType
Definition FEBasisOperations.h:1009
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:1015
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_inverseStiffnessVectorBasisType
Definition FEBasisOperations.h:1004
dftfe::uInt d_cellsBlockSize
Definition FEBasisOperations.h:1040
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessVectorBasisType
Definition FEBasisOperations.h:992
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_sqrtMassVectorCoeffType
Definition FEBasisOperations.h:984
dftfe::uInt d_nVectors
Definition FEBasisOperations.h:1038
dftfe::uInt d_nOMPThreads
Definition FEBasisOperations.h:876
std::vector< dftfe::uInt > d_nQuadsPerCell
Definition FEBasisOperations.h:1036
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:943
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:982
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:891
void clearScratchMultiVectors() const
Clears scratch multivectors.
std::vector< UpdateFlags > d_updateFlags
Definition FEBasisOperations.h:1046
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueTypeBasisData >::type, memorySpace > d_inverseMassVectorBasisTypeSinglePrec
Definition FEBasisOperations.h:974
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:628
dftfe::uInt d_quadratureID
Definition FEBasisOperations.h:1034
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > > d_JxWData
Definition FEBasisOperations.h:900
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:1023
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:1031
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellStiffnessMatrixBasisType
Definition FEBasisOperations.h:937
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > d_flattenedCellDofIndexToProcessDofIndexMap
Definition FEBasisOperations.h:890
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:1013
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellInverseSqrtMassVectorBasisType
Definition FEBasisOperations.h:954
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:909
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:970
const auto & cellSqrtMassVector() const
Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:556
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:931
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:989
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:939
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_stiffnessVectorCoeffType
Definition FEBasisOperations.h:1021
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:956
std::vector< const dealii::AffineConstraints< ValueTypeBasisData > * > * d_constraintsVector
Definition FEBasisOperations.h:878
const auto & inverseMassVector() const
Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff.
Definition FEBasisOperations.h:610
dftfe::uInt d_localSize
Definition FEBasisOperations.h:1042
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:962
void interpolate(distributedCPUVec< double > &nodalField, const dftfe::uInt dofHandlerId, const dftfe::uInt quadratureId, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureValueData, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureGradValueData, dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &quadratureHessianValueData, const bool isEvaluateGradData=false, const bool isEvaluateHessianData=false) const
void initializeShapeFunctionAndJacobianBasisData()
Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > d_cellInverseMassVectorCoeffType
Definition FEBasisOperations.h:952
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:1002
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:455
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:1045
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempCellNodalData
Definition FEBasisOperations.h:88
dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > d_cellMassVectorBasisType
Definition FEBasisOperations.h:958
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:522
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:881
dftfe::uInt d_nDofsPerCell
Definition FEBasisOperations.h:1041
std::vector< dftUtils::constraintMatrixInfo< memorySpace > > d_constraintInfo
Definition FEBasisOperations.h:875
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:928
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:960
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_JxWBasisData
Definition FEBasisOperations.h:922
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:941
void distribute(dftfe::linearAlgebra::MultiVector< ValueTypeBasisCoeff, memorySpace > &multiVector, dftfe::uInt constraintIndex=std::numeric_limits< dftfe::uInt >::max()) const
Apply constraints on given multivector.
const dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > & stiffnessVectorBasisData() const
diagonal stiffness matrix in ValueTypeBasisData
const dftfe::utils::MemoryStorage< ValueTypeBasisData, dftfe::utils::MemorySpace::HOST > & cellCentroids() const
quad point coordinates for each cell.
dftfe::utils::MemoryStorage< ValueTypeBasisCoeff, memorySpace > tempQuadratureGradientsDataNonAffine
Definition FEBasisOperations.h:89
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_shapeFunctionGradientBasisDataTranspose
Definition FEBasisOperations.h:934
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:497
const auto & shapeFunctionGradientBasisData(bool transpose=false) const
Shape function gradient values at quadrature points in ValueTypeBasisData.
Definition FEBasisOperations.h:405
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:966
std::map< dftfe::uInt, dftfe::utils::MemoryStorage< ValueTypeBasisData, memorySpace > > d_inverseJacobianBasisData
Definition FEBasisOperations.h:919
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
dealii::LinearAlgebra::distributed::Vector< elem_type, dealii::MemorySpace::Host > distributedCPUVec
Definition headers.h:92
std::uint32_t uInt
Definition TypeConfig.h:10
T type
Definition dftfeDataTypes.h:111