5 #include <amdis/GridFunctionOperator.hpp> 6 #include <amdis/Output.hpp> 7 #include <amdis/common/StaticSize.hpp> 23 template <
class LC,
class Gr
idFct>
30 static_assert( static_size_v<typename GridFct::Range> == 1,
"Expression must be of scalar type." );
37 template <
class CG,
class RN,
class CN,
class Mat>
38 void getElementMatrix(CG
const& contextGeo, RN
const& rowNode, CN
const& colNode, Mat& elementMatrix)
40 static_assert(RN::isLeaf && CN::isPower,
41 "row-node must be Leaf-Node and col-node must be a Power-Node.");
43 static const std::size_t CHILDREN = CN::CHILDREN;
44 static_assert( CHILDREN == CG::dow,
"" );
46 auto const& quad = this->getQuadratureRule(contextGeo.type(), rowNode, colNode);
47 std::size_t rowSize = rowNode.size();
48 std::size_t colSize = colNode.child(0).size();
50 using RangeFieldType =
typename CN::ChildType::LocalBasis::Traits::RangeFieldType;
51 using WorldVector = FieldVector<RangeFieldType,CG::dow>;
52 std::vector<WorldVector> colGradients;
54 for (
auto const& qp : quad) {
56 auto&& local = contextGeo.local(qp.position());
59 const auto jacobian = contextGeo.geometry().jacobianInverseTransposed(local);
62 const auto factor = Super::coefficient(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
65 auto const& shapeValues = rowNode.localBasisValuesAt(local);
68 auto const& shapeGradients = colNode.child(0).localBasisJacobiansAt(local);
71 colGradients.resize(shapeGradients.size());
73 for (std::size_t i = 0; i < colGradients.size(); ++i)
74 jacobian.mv(shapeGradients[i][0], colGradients[i]);
76 for (std::size_t i = 0; i < rowSize; ++i) {
77 const auto local_i = rowNode.localIndex(i);
78 const auto value = factor * shapeValues[i];
79 for (std::size_t j = 0; j < colSize; ++j) {
80 for (std::size_t k = 0; k < CHILDREN; ++k) {
81 const auto local_kj = colNode.child(k).localIndex(j);
82 elementMatrix[local_i][local_kj] += value * colGradients[j][k];
The base-template for GridFunctionOperators.
Definition: GridFunctionOperator.hpp:242
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
Definition: FirstOrderTestDivTrialvec.hpp:18
The main implementation of an CRTP-base class for operators using a grid-function coefficient to be u...
Definition: GridFunctionOperator.hpp:39