31 template <
class CG,
class RN,
class CN,
class Quad,
class LocalFct,
class Mat>
32 void assemble(CG
const& contextGeo, RN
const& rowNode, CN
const& colNode,
33 Quad
const& quad, LocalFct
const& localFct, Mat& elementMatrix)
const
35 static_assert(static_size_v<typename LocalFct::Range> == 1,
36 "Expression must be of scalar type.");
37 static_assert(Dune::TypeTree::Concept::UniformInnerTreeNode<RN> && Dune::TypeTree::Concept::UniformInnerTreeNode<CN>,
38 "Operator can be applied to Power-Nodes only.");
40 const bool sameFE = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
41 const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
43 if (sameFE && sameNode)
44 getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, localFct, elementMatrix);
46 getElementMatrixStandard(contextGeo, quad, rowNode, colNode, localFct, elementMatrix);
50 template <
class CG,
class QR,
class RN,
class CN,
class LocalFct,
class Mat>
51 void getElementMatrixStandard(CG
const& contextGeo, QR
const& quad,
52 RN
const& rowNode, CN
const& colNode,
53 LocalFct
const& localFct, Mat& elementMatrix)
const
55 assert( rowNode.degree() == colNode.degree() );
57 std::size_t rowSize = rowNode.child(0).size();
58 std::size_t colSize = colNode.child(0).size();
60 using RowFieldType =
typename Dune::TypeTree::Child<RN,0>::LocalBasis::Traits::RangeFieldType;
61 using RowWorldVector = FieldVector<RowFieldType,CG::dow>;
62 std::vector<RowWorldVector> rowGradients;
64 using ColFieldType =
typename Dune::TypeTree::Child<CN,0>::LocalBasis::Traits::RangeFieldType;
65 using ColWorldVector = FieldVector<ColFieldType,CG::dow>;
66 std::vector<ColWorldVector> colGradients;
68 for (
auto const& qp : quad) {
70 decltype(
auto) local = contextGeo.coordinateInElement(qp.position());
73 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
76 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
79 auto const& rowShapeGradients = rowNode.child(0).localBasisJacobiansAt(local);
80 auto const& colShapeGradients = colNode.child(0).localBasisJacobiansAt(local);
83 rowGradients.resize(rowShapeGradients.size());
85 for (std::size_t i = 0; i < rowGradients.size(); ++i)
86 jacobian.mv(rowShapeGradients[i][0], rowGradients[i]);
88 colGradients.resize(colShapeGradients.size());
90 for (std::size_t i = 0; i < colGradients.size(); ++i)
91 jacobian.mv(colShapeGradients[i][0], colGradients[i]);
93 for (std::size_t i = 0; i < rowSize; ++i) {
94 for (std::size_t j = 0; j < colSize; ++j) {
95 for (std::size_t k = 0; k < CG::dow; ++k) {
96 for (std::size_t l = 0; l < CG::dow; ++l) {
97 const auto local_ki = rowNode.child(k).localIndex(i);
98 const auto local_kj = colNode.child(l).localIndex(j);
99 elementMatrix[local_ki][local_kj] += factor * rowGradients[i][k] * colGradients[j][l];
108 template <
class CG,
class QR,
class Node,
class LocalFct,
class Mat>
109 void getElementMatrixOptimized(CG
const& contextGeo, QR
const& quad,
110 Node
const& node, Node
const& ,
111 LocalFct
const& localFct, Mat& elementMatrix)
const
113 std::size_t size = node.child(0).size();
115 using RangeFieldType =
typename Dune::TypeTree::Child<Node,0>::LocalBasis::Traits::RangeFieldType;
116 using WorldVector = FieldVector<RangeFieldType,CG::dow>;
117 std::vector<WorldVector> gradients;
119 for (
auto const& qp : quad) {
121 auto&& local = contextGeo.coordinateInElement(qp.position());
124 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
127 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
130 auto const& shapeGradients = node.child(0).localBasisJacobiansAt(local);
133 gradients.resize(shapeGradients.size());
135 for (std::size_t i = 0; i < gradients.size(); ++i)
136 jacobian.mv(shapeGradients[i][0], gradients[i]);
138 for (std::size_t k = 0; k < CG::dow; ++k) {
139 for (std::size_t l = 0; l < CG::dow; ++l) {
140 for (std::size_t i = 0; i < size; ++i) {
141 const auto local_ki = node.child(k).localIndex(i);
142 const auto value = factor * gradients[i][k];
143 elementMatrix[local_ki][local_ki] += value * gradients[i][k];
145 for (std::size_t j = i + 1; j < size; ++j) {
146 const auto local_kj = node.child(l).localIndex(j);
147 elementMatrix[local_ki][local_kj] += value * gradients[j][l];
148 elementMatrix[local_kj][local_ki] += value * gradients[j][l];