5 #include <amdis/Assembler.hpp> 6 #include <amdis/LocalOperator.hpp> 7 #include <amdis/operations/Assigner.hpp> 8 #include <amdis/typetree/Traversal.hpp> 9 #include <amdis/typetree/TreePath.hpp> 10 #include <amdis/utility/AssembleOperators.hpp> 14 template <
class GB,
class T,
class Traits>
15 template <
class ContextTag,
class Expr,
class TreePath>
17 addOperator(ContextTag contextTag, Expr
const& expr, TreePath path)
19 static_assert( Concepts::PreTreePath<TreePath>,
20 "path must be a valid treepath, or an integer/index-constant");
22 auto i = makeTreePath(path);
23 auto node = child(this->basis().localView().treeCache(), i);
25 using LocalContext =
typename ContextTag::type;
26 using Tr = DefaultAssemblerTraits<LocalContext, ElementVector>;
27 auto op = makeLocalOperator<LocalContext>(expr, this->basis().gridView());
28 auto localAssembler =
makeUniquePtr(makeAssembler<Tr>(std::move(op), node));
30 operators_[i].push(contextTag, std::move(localAssembler));
34 template <
class GB,
class T,
class Traits>
35 template <
class LocalView, REQUIRES_(Concepts::LocalView<LocalView>)>
39 elementVector_.resize(localView.
size());
42 auto const& gv = this->basis().gridView();
43 auto const& element = localView.
element();
44 auto geometry = element.geometry();
46 Traversal::forEachNode(localView.
treeCache(), [&](
auto const& node,
auto tp) {
47 auto& rhsOp = operators_[tp];
49 rhsOp.bind(element, geometry);
50 assembleOperators(gv, element, rhsOp, makeVectorAssembler(node, elementVector_));
59 template <
class GB,
class T,
class Traits>
63 auto localView = this->basis().localView();
65 this->
init(sizeInfo(*this->basis()),
true);
66 for (
auto const& element : elements(this->basis().gridView(),
typename Traits::PartitionSet{})) {
67 localView.bind(element);
68 this->assemble(localView);
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:144
Definition: Assigner.hpp:16
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:20
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:138
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:116
auto makeUniquePtr(Obj &&obj)
Create a unique_ptr by copy/move construction.
Definition: TypeTraits.hpp:107
void init(int &argc, char **&argv, std::string const &initFileName="")
Initialized the Environment for MPI.
Definition: AMDiS.hpp:29