-
Notifications
You must be signed in to change notification settings - Fork 34
/
FEBasisOperations.h
1135 lines (1005 loc) · 41.2 KB
/
FEBasisOperations.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// ---------------------------------------------------------------------
//
// Copyright (c) 2017-2022 The Regents of the University of Michigan and DFT-FE
// authors.
//
// This file is part of the DFT-FE code.
//
// The DFT-FE code is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the DFT-FE distribution.
//
// ---------------------------------------------------------------------
//
#ifndef dftfeFEBasisOperations_h
#define dftfeFEBasisOperations_h
#include <MultiVector.h>
#include <headers.h>
#include <constraintMatrixInfo.h>
#include <DeviceTypeConfig.h>
#include <BLASWrapper.h>
namespace dftfe
{
namespace basis
{
enum UpdateFlags
{
update_default = 0,
update_values = 0x0001,
update_gradients = 0x0002,
update_transpose = 0x0004,
update_quadpoints = 0x0008,
update_inversejacobians = 0x0010,
update_jxw = 0x0020,
};
inline UpdateFlags
operator|(const UpdateFlags f1, const UpdateFlags f2)
{
return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) |
static_cast<unsigned int>(f2));
}
inline UpdateFlags &
operator|=(UpdateFlags &f1, const UpdateFlags f2)
{
f1 = f1 | f2;
return f1;
}
inline UpdateFlags operator&(const UpdateFlags f1, const UpdateFlags f2)
{
return static_cast<UpdateFlags>(static_cast<unsigned int>(f1) &
static_cast<unsigned int>(f2));
}
inline UpdateFlags &
operator&=(UpdateFlags &f1, const UpdateFlags f2)
{
f1 = f1 & f2;
return f1;
}
template <typename ValueTypeBasisCoeff,
typename ValueTypeBasisData,
dftfe::utils::MemorySpace memorySpace>
class FEBasisOperations
{
protected:
mutable dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
tempCellNodalData, tempQuadratureGradientsData,
tempQuadratureGradientsDataNonAffine;
mutable dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
tempCellGradientsBlock, tempCellGradientsBlock2, tempCellValuesBlock,
tempCellMatrixBlock;
mutable dftfe::utils::MemoryStorage<dftfe::global_size_type, memorySpace>
zeroIndexVec;
mutable dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
tempCellValuesBlockCoeff;
std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
d_BLASWrapperPtr;
public:
/**
* @brief Constructor
*/
FEBasisOperations(
std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
BLASWrapperPtr);
/**
* @brief Default Destructor
*/
~FEBasisOperations() = default;
/**
* @brief Clears the FEBasisOperations internal storage.
*/
void
clear();
/**
* @brief fills required data structures for the given dofHandlerID
* @param[in] matrixFreeData MatrixFree object.
* @param[in] constraintsVector std::vector of AffineConstraints, should
* be the same vector which was passed for the construction of the given
* MatrixFree object.
* @param[in] dofHandlerID dofHandler index to be used for getting data
* from the MatrixFree object.
* @param[in] quadratureID std::vector of quadratureIDs to be used, should
* be the same IDs which were used during the construction of the given
* MatrixFree object.
*/
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);
/**
* @brief fills required data structures from another FEBasisOperations object
* @param[in] basisOperationsSrc Source FEBasisOperations object.
*/
template <dftfe::utils::MemorySpace memorySpaceSrc>
void
init(const FEBasisOperations<ValueTypeBasisCoeff,
ValueTypeBasisData,
memorySpaceSrc> &basisOperationsSrc);
/**
* @brief sets internal variables and optionally resizes internal temp storage for interpolation operations
* @param[in] vecBlockSize block size to used for operations on vectors,
* this has to be set to the exact value before any such operations are
* called.
* @param[in] cellBlockSize block size to used for cells, this has to be
* set to a value greater than or equal to the required value before any
* such operations are called
* @param[in] quadratureID Quadrature index to be used.
* @param[in] isResizeTempStorage whether to resize internal tempstorage.
*/
void
reinit(const unsigned int &vecBlockSize,
const unsigned int &cellBlockSize,
const unsigned int &quadratureID,
const bool isResizeTempStorageForInerpolation = true,
const bool isResizeTempStorageForCellMatrices = false);
dftfe::utils::MemoryStorage<dftfe::global_size_type,
dftfe::utils::MemorySpace::HOST> &
getFlattenedMapsHost();
// private:
/**
* @brief Initializes indexset maps from process level indices to cell level indices for a single vector, also initializes cell index to cellid map.
*/
void
initializeIndexMaps();
/**
* @brief Initializes indexset maps from process level indices to cell level indices for multivectors.
*/
void
initializeFlattenedIndexMaps();
/**
* @brief Initializes the constraintMatrixInfo object.
*/
void
initializeConstraints();
/**
* @brief Reinitializes the constraintMatrixInfo object.
*/
void
reinitializeConstraints(
std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
&constraintsVector);
/**
* @brief Constructs the MPIPatternP2P object.
*/
void
initializeMPIPattern();
/**
* @brief Fill the shape function data and jacobian data in the ValueTypeBasisCoeff datatype.
*/
void
initializeShapeFunctionAndJacobianData();
/**
* @brief Fill the shape function data and jacobian data in the ValueTypeBasisData datatype.
*/
void
initializeShapeFunctionAndJacobianBasisData();
/**
* @brief Computes the cell-level stiffness matrix.
*/
void
computeCellStiffnessMatrix(const unsigned int quadratureID,
const unsigned int cellsBlockSize,
const bool basisType = false,
const bool ceoffType = true);
void
computeCellMassMatrix(const unsigned int quadratureID,
const unsigned int cellsBlockSize,
const bool basisType = false,
const bool ceoffType = true);
void
computeWeightedCellMassMatrix(
const std::pair<unsigned int, unsigned int> cellRangeTotal,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &weights,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
&weightedCellMassMatrix) const;
void
computeWeightedCellNjGradNiMatrix(
const std::pair<unsigned int, unsigned int> cellRangeTotal,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &weights,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
&weightedCellNjGradNiMatrix) const;
void
computeWeightedCellNjGradNiPlusNiGradNjMatrix(
const std::pair<unsigned int, unsigned int> cellRangeTotal,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &weights,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
&weightedCellNjGradNiPlusNiGradNjMatrix) const;
void
computeInverseSqrtMassVector(const bool basisType = true,
const bool ceoffType = false);
/**
* @brief Computes the stiffness matrix
* \grad Ni . \grad Ni.
*/
void
computeStiffnessVector(const bool basisType = true,
const bool ceoffType = false);
/**
* @brief Resizes the internal temp storage to be sufficient for the vector and cell block sizes provided in reinit.
*/
void
resizeTempStorage(const bool isResizeTempStorageForInerpolation,
const bool isResizeTempStorageForCellMatrices);
/**
* @brief Number of quadrature points per cell for the quadratureID set in reinit.
*/
unsigned int
nQuadsPerCell() const;
/**
* @brief Number of DoFs per cell for the dofHandlerID set in init.
*/
unsigned int
nDofsPerCell() const;
/**
* @brief Number of locally owned cells on the current processor.
*/
unsigned int
nCells() const;
/**
* @brief Number of DoFs on the current processor, locally owned + ghosts.
*/
unsigned int
nRelaventDofs() const;
/**
* @brief Number of locally owned DoFs on the current processor.
*/
unsigned int
nOwnedDofs() const;
/**
* @brief Shape function values at quadrature points.
* @param[in] transpose if false the the data is indexed as [iQuad *
* d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
* d_nQuadsPerCell + iQuad].
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> &
shapeFunctionData(bool transpose = false) const;
/**
* @brief Shape function gradient values at quadrature points.
* @param[in] transpose if false the the data is indexed as [iDim *
* d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
* if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
* iNode * d_nQuadsPerCell + iQuad].
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> &
shapeFunctionGradientData(bool transpose = false) const;
/**
* @brief Inverse Jacobian matrices, for cartesian cells returns the
* diagonal elements of the inverse Jacobian matrices for each cell, for
* affine cells returns the 3x3 inverse Jacobians for each cell otherwise
* returns the 3x3 inverse Jacobians at each quad point for each cell.
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> &
inverseJacobians() const;
/**
* @brief determinant of Jacobian times the quadrature weight at each
* quad point for each cell.
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace> &
JxW() const;
/**
* @brief quad point coordinates for each cell.
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData,
dftfe::utils::MemorySpace::HOST> &
quadPoints() const;
/**
* @brief Shape function values at quadrature points in ValueTypeBasisData.
* @param[in] transpose if false the the data is indexed as [iQuad *
* d_nDofsPerCell + iNode] and if true it is indexed as [iNode *
* d_nQuadsPerCell + iQuad].
*/
const auto &
shapeFunctionBasisData(bool transpose = false) const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return transpose ?
d_shapeFunctionDataTranspose.find(d_quadratureID)->second :
d_shapeFunctionData.find(d_quadratureID)->second;
}
else
{
return transpose ?
d_shapeFunctionBasisDataTranspose.find(d_quadratureID)
->second :
d_shapeFunctionBasisData.find(d_quadratureID)->second;
}
}
/**
* @brief Shape function gradient values at quadrature points in ValueTypeBasisData.
* @param[in] transpose if false the the data is indexed as [iDim *
* d_nQuadsPerCell * d_nDofsPerCell + iQuad * d_nDofsPerCell + iNode] and
* if true it is indexed as [iDim * d_nQuadsPerCell * d_nDofsPerCell +
* iNode * d_nQuadsPerCell + iQuad].
*/
const auto &
shapeFunctionGradientBasisData(bool transpose = false) const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return transpose ?
d_shapeFunctionGradientDataTranspose.find(d_quadratureID)
->second :
d_shapeFunctionGradientData.find(d_quadratureID)->second;
}
else
{
return transpose ?
d_shapeFunctionGradientBasisDataTranspose
.find(d_quadratureID)
->second :
d_shapeFunctionGradientBasisData.find(d_quadratureID)
->second;
}
}
/**
* @brief Inverse Jacobian matrices in ValueTypeBasisData, for cartesian cells returns the
* diagonal elements of the inverse Jacobian matrices for each cell, for
* affine cells returns the 3x3 inverse Jacobians for each cell otherwise
* returns the 3x3 inverse Jacobians at each quad point for each cell.
*/
const auto &
inverseJacobiansBasisData() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_inverseJacobianData
.find(areAllCellsAffine ? 0 : d_quadratureID)
->second;
}
else
{
return d_inverseJacobianBasisData
.find(areAllCellsAffine ? 0 : d_quadratureID)
->second;
}
}
/**
* @brief determinant of Jacobian times the quadrature weight in ValueTypeBasisData at each
* quad point for each cell.
*/
const auto &
JxWBasisData() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_JxWData.find(d_quadratureID)->second;
}
else
{
return d_JxWBasisData.find(d_quadratureID)->second;
}
}
/**
* @brief Cell level stiffness matrix in ValueTypeBasisCoeff
*/
const auto &
cellStiffnessMatrix() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_cellStiffnessMatrixBasisType;
}
else
{
return d_cellStiffnessMatrixCoeffType;
}
}
/**
* @brief Cell level stiffness matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
cellStiffnessMatrixBasisData() const;
/**
* @brief Cell level mass matrix in ValueTypeBasisCoeff
*/
const auto &
cellMassMatrix() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_cellMassMatrixBasisType;
}
else
{
return d_cellMassMatrixCoeffType;
}
}
/**
* @brief Cell level mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
cellMassMatrixBasisData() const;
/**
* @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
cellInverseSqrtMassVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_cellInverseSqrtMassVectorBasisType;
}
else
{
return d_cellInverseSqrtMassVectorCoeffType;
}
}
/**
* @brief Cell level inverse diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
cellInverseMassVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_cellInverseMassVectorBasisType;
}
else
{
return d_cellInverseMassVectorCoeffType;
}
}
/**
* @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
cellSqrtMassVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_cellSqrtMassVectorBasisType;
}
else
{
return d_cellSqrtMassVectorCoeffType;
}
}
/**
* @brief Cell level diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
cellMassVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_cellMassVectorBasisType;
}
else
{
return d_cellMassVectorCoeffType;
}
}
/**
* @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
inverseSqrtMassVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_inverseSqrtMassVectorBasisType;
}
else
{
return d_inverseSqrtMassVectorCoeffType;
}
}
/**
* @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
inverseMassVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_inverseMassVectorBasisType;
}
else
{
return d_inverseMassVectorCoeffType;
}
}
/**
* @brief sqrt diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
sqrtMassVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_sqrtMassVectorBasisType;
}
else
{
return d_sqrtMassVectorCoeffType;
}
}
/**
* @brief diagonal mass matrix in ValueTypeBasisCoeff
*/
const auto &
massVector() const
{
if constexpr (std::is_same<ValueTypeBasisCoeff,
ValueTypeBasisData>::value)
{
return d_massVectorBasisType;
}
else
{
return d_massVectorCoeffType;
}
}
/**
* @brief Cell level inverse sqrt diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
cellInverseSqrtMassVectorBasisData() const;
/**
* @brief Cell level inverse diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
cellInverseMassVectorBasisData() const;
/**
* @brief Cell level sqrt diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
cellSqrtMassVectorBasisData() const;
/**
* @brief Cell level diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
cellMassVectorBasisData() const;
/**
* @brief Inverse sqrt diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
inverseSqrtMassVectorBasisData() const;
/**
* @brief Inverse diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
inverseMassVectorBasisData() const;
/**
* @brief sqrt diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
sqrtMassVectorBasisData() const;
/**
* @brief diagonal mass matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
massVectorBasisData() const;
/**
* @brief diagonal stiffness matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
stiffnessVectorBasisData() const;
/**
* @brief diagonal inverse stiffness matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
inverseStiffnessVectorBasisData() const;
/**
* @brief diagonal inverse sqrt stiffness matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
inverseSqrtStiffnessVectorBasisData() const;
/**
* @brief diagonal sqrt stiffness matrix in ValueTypeBasisData
*/
const dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace> &
sqrtStiffnessVectorBasisData() const;
/**
* @brief returns 2 if all cells on current processor are Cartesian,
* 1 if all cells on current processor are affine and 0 otherwise.
*/
unsigned int
cellsTypeFlag() const;
/**
* @brief returns the deal.ii cellID corresponing to given cell Index.
* @param[in] iElem cell Index
*/
dealii::CellId
cellID(const unsigned int iElem) const;
/**
* @brief returns the cell index corresponding to given deal.ii cellID.
* @param[in] iElem cell Index
*/
unsigned int
cellIndex(const dealii::CellId cellid) const;
/**
* @brief Creates a multivector.
* @param[in] blocksize Number of vectors in the multivector.
* @param[out] multiVector the created multivector.
*/
void
createMultiVector(
const unsigned int blocksize,
dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff, memorySpace>
&multiVector) const;
/**
* @brief Creates scratch multivectors.
* @param[in] vecBlockSize Number of vectors in the multivector.
* @param[out] numMultiVecs number of scratch multivectors needed with
* this vecBlockSize.
*/
void
createScratchMultiVectors(const unsigned int vecBlockSize,
const unsigned int numMultiVecs = 1) const;
/**
* @brief Creates single precision scratch multivectors.
* @param[in] vecBlockSize Number of vectors in the multivector.
* @param[out] numMultiVecs number of scratch multivectors needed with
* this vecBlockSize.
*/
void
createScratchMultiVectorsSinglePrec(
const unsigned int vecBlockSize,
const unsigned int numMultiVecs = 1) const;
/**
* @brief Clears scratch multivectors.
*/
void
clearScratchMultiVectors() const;
/**
* @brief Gets scratch multivectors.
* @param[in] vecBlockSize Number of vectors in the multivector.
* @param[out] numMultiVecs index of the multivector among those with the
* same vecBlockSize.
*/
dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff, memorySpace> &
getMultiVector(const unsigned int vecBlockSize,
const unsigned int index = 0) const;
/**
* @brief Gets single precision scratch multivectors.
* @param[in] vecBlockSize Number of vectors in the multivector.
* @param[out] numMultiVecs index of the multivector among those with the
* same vecBlockSize.
*/
dftfe::linearAlgebra::MultiVector<
typename dftfe::dataTypes::singlePrecType<ValueTypeBasisCoeff>::type,
memorySpace> &
getMultiVectorSinglePrec(const unsigned int vecBlockSize,
const unsigned int index = 0) const;
/**
* @brief Apply constraints on given multivector.
* @param[inout] multiVector the given multivector.
*/
void
distribute(dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
memorySpace> &multiVector,
unsigned int constraintIndex =
std::numeric_limits<unsigned int>::max()) const;
/**
* @brief Return the underlying deal.II matrixfree object.
*/
const dealii::MatrixFree<3, ValueTypeBasisData> &
matrixFreeData() const;
/**
* @brief Return the underlying deal.II dofhandler object.
*/
const dealii::DoFHandler<3> &
getDofHandler() const;
std::vector<dftUtils::constraintMatrixInfo<memorySpace>> d_constraintInfo;
unsigned int d_nOMPThreads;
std::vector<const dealii::AffineConstraints<ValueTypeBasisData> *>
* d_constraintsVector;
const dealii::MatrixFree<3, ValueTypeBasisData> *d_matrixFreeDataPtr;
dftfe::utils::MemoryStorage<dftfe::global_size_type,
dftfe::utils::MemorySpace::HOST>
d_cellDofIndexToProcessDofIndexMap;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisData,
dftfe::utils::MemorySpace::HOST>>
d_quadPoints;
dftfe::utils::MemoryStorage<dftfe::global_size_type, memorySpace>
d_flattenedCellDofIndexToProcessDofIndexMap;
std::vector<dealii::CellId> d_cellIndexToCellIdMap;
std::map<dealii::CellId, unsigned int> d_cellIdToCellIndexMap;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>>
d_inverseJacobianData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>>
d_JxWData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>>
d_shapeFunctionData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>>
d_shapeFunctionGradientDataInternalLayout;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>>
d_shapeFunctionGradientData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>>
d_shapeFunctionDataTranspose;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>>
d_shapeFunctionGradientDataTranspose;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>>
d_inverseJacobianBasisData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>>
d_JxWBasisData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>>
d_shapeFunctionBasisData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>>
d_shapeFunctionGradientBasisData;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>>
d_shapeFunctionBasisDataTranspose;
std::map<unsigned int,
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>>
d_shapeFunctionGradientBasisDataTranspose;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellStiffnessMatrixBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellStiffnessMatrixCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellMassMatrixBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellMassMatrixCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellInverseMassVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellInverseMassVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellInverseSqrtMassVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellInverseSqrtMassVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellMassVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellMassVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellSqrtMassVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellSqrtMassVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_massVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_massVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_inverseMassVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_inverseMassVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_inverseSqrtMassVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_inverseSqrtMassVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_sqrtMassVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_sqrtMassVectorCoeffType;
mutable std::map<
unsigned int,
std::vector<
dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff, memorySpace>>>
scratchMultiVectors;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellStiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellInverseStiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellSqrtStiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_cellInverseSqrtStiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_inverseSqrtStiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_sqrtStiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_inverseStiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisData, memorySpace>
d_stiffnessVectorBasisType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellInverseStiffnessVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellInverseSqrtStiffnessVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellStiffnessVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_cellSqrtStiffnessVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_inverseSqrtStiffnessVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_sqrtStiffnessVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_stiffnessVectorCoeffType;
dftfe::utils::MemoryStorage<ValueTypeBasisCoeff, memorySpace>
d_inverseStiffnessVectorCoeffType;
mutable std::map<
unsigned int,
std::vector<dftfe::linearAlgebra::MultiVector<
typename dftfe::dataTypes::singlePrecType<ValueTypeBasisCoeff>::type,
memorySpace>>>
scratchMultiVectorsSinglePrec;
std::vector<unsigned int> d_quadratureIDsVector;
unsigned int d_quadratureID;
unsigned int d_quadratureIndex;
std::vector<unsigned int> d_nQuadsPerCell;
unsigned int d_dofHandlerID;
unsigned int d_nVectors;
unsigned int d_nCells;
unsigned int d_cellsBlockSize;
unsigned int d_nDofsPerCell;
unsigned int d_localSize;
unsigned int d_locallyOwnedSize;
bool areAllCellsAffine;
bool areAllCellsCartesian;
std::vector<UpdateFlags> d_updateFlags;
std::shared_ptr<const utils::mpi::MPIPatternP2P<memorySpace>>
mpiPatternP2P;
/**
* @brief Interpolate process level nodal data to cell level quadrature data.
* @param[in] nodalData process level nodal data, the multivector should
* already have ghost data and constraints should have been applied.
* @param[out] quadratureValues Cell level quadrature values, indexed by
* [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
* @param[out] quadratureGradients Cell level quadrature gradients,
* indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *
* d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
*/
void
interpolate(dftfe::linearAlgebra::MultiVector<ValueTypeBasisCoeff,
memorySpace> &nodalData,
ValueTypeBasisCoeff *quadratureValues,
ValueTypeBasisCoeff *quadratureGradients = NULL) const;
// FIXME Untested function
/**
* @brief Integrate cell level quadrature data times shape functions to process level nodal data.
* @param[in] quadratureValues Cell level quadrature values, indexed by
* [iCell * d_nQuadsPerCell * d_nVectors + iQuad * d_nVectors + iVec].
* @param[in] quadratureGradients Cell level quadrature gradients,
* indexed by [iCell * 3 * d_nQuadsPerCell * d_nVectors + iDim *