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 using expr_value_type =
typename LocalFct::Range;
36 static_assert(static_size_v<expr_value_type> == 1 ||
37 (static_num_rows_v<expr_value_type> == CG::dow &&
38 static_num_cols_v<expr_value_type> == CG::dow),
39 "Expression must be of scalar or matrix type." );
40 static_assert(Dune::TypeTree::Concept::UniformInnerTreeNode<RN> && Dune::TypeTree::Concept::UniformInnerTreeNode<CN>,
41 "Operator can be applied to Power-Nodes only.");
42 assert(rowNode.degree() == colNode.degree());
44 const bool sameFE = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
45 const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
47 using Category = ValueCategory_t<expr_value_type>;
49 if (sameFE && sameNode && std::is_same_v<Category,tag::scalar>)
50 getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, localFct, elementMatrix, Category{});
52 getElementMatrixStandard(contextGeo, quad, rowNode, colNode, localFct, elementMatrix, Category{});
57 template <
class CG,
class QR,
class RN,
class CN,
class LocalFct,
class Mat>
58 void getElementMatrixStandard(CG
const& contextGeo, QR
const& quad,
59 RN
const& rowNode, CN
const& colNode,
60 LocalFct
const& localFct, Mat& elementMatrix,
tag::scalar)
const
62 std::size_t rowSize = rowNode.child(0).size();
63 std::size_t colSize = colNode.child(0).size();
65 for (
auto const& qp : quad) {
67 auto&& local = contextGeo.coordinateInElement(qp.position());
70 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
72 auto const& rowShapeValues = rowNode.child(0).localBasisValuesAt(local);
73 auto const& colShapeValues = colNode.child(0).localBasisValuesAt(local);
75 for (std::size_t i = 0; i < rowSize; ++i) {
76 const auto value = factor * rowShapeValues[i];
78 for (std::size_t j = 0; j < colSize; ++j) {
79 const auto value0 = value * colShapeValues[j];
81 for (std::size_t k = 0; k < rowNode.degree(); ++k) {
82 const auto local_ki = rowNode.child(k).localIndex(i);
83 const auto local_kj = colNode.child(k).localIndex(j);
85 elementMatrix[local_ki][local_kj] += value0;
93 template <
class CG,
class QR,
class RN,
class CN,
class LocalFct,
class Mat>
94 void getElementMatrixOptimized(CG
const& contextGeo, QR
const& quad,
95 RN
const& node, CN
const& ,
96 LocalFct
const& localFct, Mat& elementMatrix,
tag::scalar)
const
98 std::size_t size = node.child(0).size();
100 for (
auto const& qp : quad) {
102 auto&& local = contextGeo.coordinateInElement(qp.position());
105 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
107 auto const& shapeValues = node.child(0).localBasisValuesAt(local);
109 for (std::size_t i = 0; i < size; ++i) {
110 const auto value = factor * shapeValues[i];
112 for (std::size_t k = 0; k < node.degree(); ++k) {
113 const auto local_ki = node.child(k).localIndex(i);
114 elementMatrix[local_ki][local_ki] += value * shapeValues[i];
117 for (std::size_t j = i+1; j < size; ++j) {
118 const auto value0 = value * shapeValues[j];
120 for (std::size_t k = 0; k < node.degree(); ++k) {
121 const auto local_ki = node.child(k).localIndex(i);
122 const auto local_kj = node.child(k).localIndex(j);
124 elementMatrix[local_ki][local_kj] += value0;
125 elementMatrix[local_kj][local_ki] += value0;
132 template <
class CG,
class QR,
class RN,
class CN,
class LocalFct,
class Mat>
133 void getElementMatrixStandard(CG
const& contextGeo, QR
const& quad,
134 RN
const& rowNode, CN
const& colNode,
135 LocalFct
const& localFct, Mat& elementMatrix,
tag::matrix)
const
137 using expr_value_type [[maybe_unused]] =
typename LocalFct::Range;
138 assert(static_num_rows_v<expr_value_type> == rowNode.degree() &&
139 static_num_cols_v<expr_value_type> == rowNode.degree());
141 std::size_t rowSize = rowNode.child(0).size();
142 std::size_t colSize = colNode.child(0).size();
144 for (
auto const& qp : quad) {
146 auto&& local = contextGeo.coordinateInElement(qp.position());
149 const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
150 const auto exprValue = localFct(local);
152 auto const& rowShapeValues = rowNode.child(0).localBasisValuesAt(local);
153 auto const& colShapeValues = colNode.child(0).localBasisValuesAt(local);
155 for (std::size_t i = 0; i < rowSize; ++i) {
156 const auto value0 = factor * rowShapeValues[i];
158 for (std::size_t j = 0; j < colSize; ++j) {
159 const auto value = value0 * colShapeValues[j];
160 const auto mat = exprValue * value;
162 for (std::size_t k0 = 0; k0 < rowNode.degree(); ++k0) {
163 const auto local_ki = rowNode.child(k0).localIndex(i);
164 for (std::size_t k1 = 0; k1 < rowNode.degree(); ++k1) {
165 const auto local_kj = colNode.child(k1).localIndex(j);
167 elementMatrix[local_ki][local_kj] += mat[k0][k1];
175 template <
class CG,
class QR,
class RN,
class CN,
class LocalFct,
class Mat>
176 void getElementMatrixOptimized(CG
const&, QR
const&, RN
const&, CN
const&,
179 DUNE_THROW(Dune::NotImplemented,
180 "Optimized LocalOperator with matrix coefficients not implemented.");