AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
ZeroOrderTestTrialvec.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7#include <amdis/common/FieldMatVec.hpp>
8
9#include <dune/common/typetree/nodeconcepts.hh>
10
11namespace AMDiS
12{
18 namespace tag
19 {
20 struct test_trialvec {};
21 }
22
23
26 {
27 public:
29
30 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
31 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
32 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
33 {
34 static_assert(static_size_v<typename LocalFct::Range> == CG::dow,
35 "Expression must be of vector type." );
36 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::UniformInnerTreeNode<CN>,
37 "RN must be Leaf-Node and CN must be a Power-Node.");
38
39 assert(colNode.degree() == CG::dow);
40
41 std::size_t rowSize = rowNode.size();
42 std::size_t colSize = colNode.child(0).size();
43
44 for (auto const& qp : quad) {
45 // Position of the current quadrature point in the reference element
46 auto&& local = contextGeo.coordinateInElement(qp.position());
47
48 // The multiplicative factor in the integral transformation formula
49 const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
50 const auto b = localFct(local);
51
52 auto const& rowShapeValues = rowNode.localBasisValuesAt(local);
53 auto const& colShapeValues = colNode.child(0).localBasisValuesAt(local);
54
55 for (std::size_t i = 0; i < rowSize; ++i) {
56 const auto local_i = rowNode.localIndex(i);
57
58 for (std::size_t j = 0; j < colSize; ++j) {
59 const auto value = b * (factor * rowShapeValues[i] * colShapeValues[j]);
60
61 for (std::size_t k = 0; k < colNode.degree(); ++k) {
62 const auto local_kj = colNode.child(k).localIndex(j);
63 elementMatrix[local_i][local_kj] += Dune::at(value,k);
64 }
65 }
66 }
67 }
68 }
69 };
70
71 template <class LC>
72 struct GridFunctionOperatorRegistry<tag::test_trialvec, LC>
73 {
74 static constexpr int degree = 0;
76 };
77
80} // end namespace AMDiS
zero-order operator
Definition ZeroOrderTestTrialvec.hpp:26
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition ZeroOrderTestTrialvec.hpp:20