5 #include <amdis/GridFunctionOperator.hpp> 6 #include <amdis/common/StaticSize.hpp> 7 #include <amdis/utility/LocalToGlobalAdapter.hpp> 40 template <
class CG,
class RN,
class CN,
class Quad,
class LocalFct,
class Mat>
41 void assemble(CG
const& contextGeo, RN
const& rowNode, CN
const& colNode,
42 Quad
const& quad, LocalFct
const& localFct, Mat& elementMatrix)
const 44 static_assert(static_size_v<typename LocalFct::Range> == 1,
45 "Expression must be of scalar type.");
46 static_assert(RN::isLeaf && CN::isLeaf,
47 "Operator can be applied to Leaf-Nodes only.");
49 std::size_t rowSize = rowNode.size();
50 std::size_t colSize = colNode.size();
53 for (
auto const& qp : quad) {
55 decltype(
auto) local = contextGeo.coordinateInElement(qp.position());
58 auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
59 factor *= localFct(local);
62 auto const& shapeValues = rowNode.localBasisValuesAt(local);
65 auto const& colPartial = colLocalBasis.partialsAt(local, comp_);
67 for (std::size_t j = 0; j < colSize; ++j) {
68 const auto local_j = colNode.localIndex(j);
69 const auto value = factor * colPartial[j];
70 for (std::size_t i = 0; i < rowSize; ++i) {
71 const auto local_i = rowNode.localIndex(i);
72 elementMatrix[local_i][local_j] += value * shapeValues[i];
82 static constexpr
int degree = 1;
Convert a simple (scalar) local basis into a global basis.
Definition: LocalToGlobalAdapter.hpp:64
Definition: FirstOrderTestPartialTrial.hpp:23
Definition: AdaptBase.hpp:6
Definition: FirstOrderTestPartialTrial.hpp:18
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition: GridFunctionOperator.hpp:235
auto fot(Expr &&expr, tag::partial_trial t, int quadOrder=-1)
Create a first-order term with derivative on trial-function.
Definition: FirstOrderTestPartialTrial.hpp:89
first-order operator
Definition: FirstOrderTestPartialTrial.hpp:31
Registry to specify a tag for each implementation type.
Definition: GridFunctionOperator.hpp:216