AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
FirstOrderPartialTest.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7#include <amdis/utility/LocalToGlobalAdapter.hpp>
8
9#include <dune/common/typetree/nodeconcepts.hh>
10
11namespace AMDiS
12{
18 namespace tag
19 {
21 {
22 std::size_t comp;
23 };
24 }
25
26
29 {
30 int comp_;
31
32 public:
34 : comp_(tag.comp)
35 {}
36
37 template <class CG, class Node, class Quad, class LocalFct, class Vec>
38 void assemble(CG const& contextGeo, Node const& node, Quad const& quad,
39 LocalFct const& localFct, Vec& elementVector) const
40 {
41 static_assert(static_size_v<typename LocalFct::Range> == 1, "Expression must be of scalar type." );
42 static_assert(Dune::TypeTree::Concept::LeafTreeNode<Node>, "Operator can be applied to Leaf-Nodes only.");
43
44 LocalToGlobalBasisAdapter localBasis(node, contextGeo.elementGeometry());
45
46 for (auto const& qp : quad) {
47 // Position of the current quadrature point in the reference element
48 decltype(auto) local = contextGeo.coordinateInElement(qp.position());
49
50 // The multiplicative factor in the integral transformation formula
51 auto dx = contextGeo.integrationElement(qp.position()) * qp.weight();
52 dx *= localFct(local);
53
54 // Compute the shape function gradients on the real element
55 auto const& partial = localBasis.partialsAt(local, comp_);
56
57 for (std::size_t j = 0; j < localBasis.size(); ++j) {
58 const auto local_j = node.localIndex(j);
59 elementVector[local_j] += dx * partial[j];
60 }
61 }
62 }
63 };
64
65 template <class LC>
66 struct GridFunctionOperatorRegistry<tag::partialtest, LC>
67 {
68 static constexpr int degree = 1;
70 };
71
74} // end namespace AMDiS
first-order operator
Definition FirstOrderPartialTest.hpp:29
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
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition FirstOrderPartialTest.hpp:21