5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
27 template <
class CG,
class RN,
class CN,
class Quad,
class LocalFct,
class Mat>
28 void assemble(CG
const& contextGeo, RN
const& rowNode, CN
const& colNode,
29 Quad
const& quad, LocalFct
const& localFct, Mat& elementMatrix)
const
31 static_assert(static_size_v<typename LocalFct::Range> == 1,
32 "Expression must be of scalar type.");
33 static_assert(RN::isPower && CN::isLeaf,
34 "col-node must be Leaf-Node and row-node must be a Power-Node.");
36 assert( rowNode.degree() == CG::dow );
38 std::size_t rowSize = rowNode.child(0).size();
39 std::size_t colSize = colNode.size();
41 using RangeFieldType =
typename CN::LocalBasis::Traits::RangeFieldType;
42 using WorldVector = FieldVector<RangeFieldType,CG::dow>;
43 std::vector<WorldVector> colGradients;
45 for (
auto const& qp : quad) {
47 decltype(
auto) local = contextGeo.coordinateInElement(qp.position());
50 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
53 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
56 auto const& shapeValues = rowNode.child(0).localBasisValuesAt(local);
59 auto const& shapeGradients = colNode.localBasisJacobiansAt(local);
62 colGradients.resize(shapeGradients.size());
64 for (std::size_t i = 0; i < colGradients.size(); ++i)
65 jacobian.mv(shapeGradients[i][0], colGradients[i]);
67 for (std::size_t i = 0; i < rowSize; ++i) {
68 const auto value = factor * shapeValues[i];
69 for (std::size_t j = 0; j < colSize; ++j) {
70 const auto local_j = colNode.localIndex(j);
71 for (std::size_t k = 0; k < rowNode.degree(); ++k) {
72 const auto local_ki = rowNode.child(k).localIndex(i);
73 elementMatrix[local_ki][local_j] += value * colGradients[j][k];
85 static constexpr int degree = 1;
first-order operator
Definition: FirstOrderTestvecGradTrial.hpp:23
Registry to specify a tag for each implementation type.
Definition: GridFunctionOperator.hpp:216
Definition: FirstOrderTestvecGradTrial.hpp:17