AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
SecondOrderPartialTestPartialTrial.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7
8#include <dune/common/typetree/nodeconcepts.hh>
9
10namespace AMDiS
11{
17 namespace tag
18 {
20 {
21 std::size_t comp_test; // i
22 std::size_t comp_trial; // j
23 };
24 }
25
26
29 {
30 int compTest_;
31 int compTrial_;
32
33 public:
35 : compTest_(tag.comp_test)
36 , compTrial_(tag.comp_trial)
37 {}
38
39 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
40 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
41 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
42 {
43 static_assert(static_size_v<typename LocalFct::Range> == 1,
44 "Expression must be of scalar type." );
45 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
46 "Operator can be applied to Leaf-Nodes only.");
47
48 LocalToGlobalBasisAdapter rowLocalBasis(rowNode, contextGeo.elementGeometry());
49 LocalToGlobalBasisAdapter colLocalBasis(colNode, contextGeo.elementGeometry());
50
51 for (auto const& qp : quad) {
52 // Position of the current quadrature point in the reference element
53 decltype(auto) local = contextGeo.coordinateInElement(qp.position());
54
55 // The multiplicative factor in the integral transformation formula
56 auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
57 factor *= localFct(local);
58
59 auto const& rowPartial = rowLocalBasis.partialsAt(local, compTest_);
60 auto const& colPartial = colLocalBasis.partialsAt(local, compTrial_);
61
62 for (std::size_t j = 0; j < colLocalBasis.size(); ++j) {
63 const auto local_j = colNode.localIndex(j);
64 const auto value = factor * colPartial[j];
65 for (std::size_t i = 0; i < rowLocalBasis.size(); ++i) {
66 const auto local_i = rowNode.localIndex(i);
67 elementMatrix[local_i][local_j] += value * rowPartial[i];
68 }
69 }
70 }
71
72 }
73 };
74
75 template <class LC>
76 struct GridFunctionOperatorRegistry<tag::partialtest_partialtrial, LC>
77 {
78 static constexpr int degree = 2;
80 };
81
82
84 template <class Expr>
85 auto sot_ij(Expr&& expr, std::size_t comp_test, std::size_t comp_trial, int quadOrder = -1)
86 {
87 return makeOperator(tag::partialtest_partialtrial{comp_test, comp_trial}, FWD(expr), quadOrder);
88 }
89
92} // end namespace AMDiS
Convert a simple (scalar) local basis into a global basis.
Definition LocalToGlobalAdapter.hpp:65
auto const & partialsAt(typename Traits::DomainLocal const &x, std::size_t comp) const
Definition LocalToGlobalAdapter.hpp:175
std::size_t size() const
Return the number of local basis functions.
Definition LocalToGlobalAdapter.hpp:96
second-order operator
Definition SecondOrderPartialTestPartialTrial.hpp:29
auto sot_ij(Expr &&expr, std::size_t comp_test, std::size_t comp_trial, int quadOrder=-1)
Create a second-order term of partial derivatives.
Definition SecondOrderPartialTestPartialTrial.hpp:85
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition GridFunctionOperator.hpp:235
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition SecondOrderPartialTestPartialTrial.hpp:20