13#include <dune/common/float_cmp.hh>
14#include <dune/grid/utility/hierarchicsearch.hh>
33 auto fT = [&]([[maybe_unused]]
int nodeNumber,
35 lagrangeNodeGlobalCoords.emplace_back(localView.element().geometry().global(localCoordinate));
60template <
int worldDim,
typename Basis, Dune::FloatCmp::CmpStyle cmpStyle = Dune::FloatCmp::CmpStyle::absolute>
62 constexpr double tol = 1e-8;
63 using LocalView = std::remove_cvref_t<
decltype(
basis.localView())>;
64 using MultiIndex =
typename LocalView::MultiIndex;
65 using Element =
typename LocalView::Element;
66 using Tree = LocalView::Tree;
67 static constexpr bool isScalar = Tree::isLeaf;
68 static constexpr bool isPower = Tree::isPower;
70 if constexpr (Tree::isComposite) {
71 Dune::Hybrid::forEach(
72 Dune::Hybrid::integralRange(Dune::index_constant<std::tuple_size_v<typename Basis::PreBasis::SubPreBases>>()),
75 not Tree::template Child<i>::Type::isComposite,
76 "globalIndexFromGlobalPosition is not implemented to handle a composite basis within a composite basis.");
80 static constexpr std::size_t numChildren = Impl::PreBasisInfo<Tree>::size;
81 constexpr int myDim = Element::mydimension;
82 Dune::HierarchicSearch hSearch(
basis.gridView().grid(),
basis.gridView().indexSet());
83 std::conditional_t<isScalar, std::optional<MultiIndex>, std::optional<std::array<MultiIndex, numChildren>>>
85 const auto& ele = hSearch.findEntity(pos);
86 auto localView =
basis.localView();
88 const auto geo = localView.element().geometry();
89 const auto& node = [&]() {
90 if constexpr (isScalar or isPower)
91 return localView.tree();
93 return localView.tree().template child<0>();
96 Impl::checkLagrangeNode<isScalar>(node);
101 if constexpr (isScalar)
102 globalIndices.value() = localView.index(node.localIndex(nodeNumber));
104 for (
int j = 0; j < numChildren; j++)
105 globalIndices.value()[j] = localView.index(node.child(j).localIndex(nodeNumber));
113 DUNE_THROW(Dune::InvalidStateException,
"No Lagrange node found at the given position in the grid.");
125template <
typename FE>
129 const auto& referenceElement = Dune::referenceElement<double, dim>(element.type());
130 const int numberOfVertices = referenceElement.size(codim);
132 auto getPosition = [=](
const int i) {
return referenceElement.position(i, codim); };
133 return std::views::transform(std::views::iota(0, numberOfVertices), getPosition);
143template <
typename FE>
172 auto& dv = assembler->dirichletValues();
173 auto newInc = Eigen::VectorXd::Zero(dv.size()).eval();
174 dv.evaluateInhomogeneousBoundaryConditionDerivative(newInc, 1.0);
176 Eigen::VectorXd F_dirichlet;
177 F_dirichlet = K * newInc;
180 assembler->dirichletValues().setZeroAtConstrainedDofs(F_dirichlet);
182 F_dirichlet = assembler->createReducedVector(F_dirichlet);
Implementation of creating a flat basis from a possibly blocked basis.
Contains functions to traverse through a tree to its different nodes.
auto referenceElementSubEntityPositions(FE &fe, int codim)
A function to obtain the local coordinates of subentities of an FiniteElement.
Definition: functionhelper.hh:126
auto globalIndexFromGlobalPosition(const Basis &basis, const Dune::FieldVector< double, worldDim > &pos)
A helper function to obtain the global index from the global positions for a Lagrange node.
Definition: functionhelper.hh:61
auto referenceElementVertexPositions(FE &fe)
A function to obtain the local coordinates the vertices of an FiniteElement.
Definition: functionhelper.hh:144
void obtainLagrangeGlobalNodePositions(const LV &localView, std::vector< Dune::FieldVector< double, LV::Element::Geometry::coorddimension > > &lagrangeNodeGlobalCoords)
A function to obtain the global positions of the nodes of an element with Lagrangian basis.
Definition: functionhelper.hh:30
void globalIndices(const FiniteElement &fe, std::vector< typename FiniteElement::LocalView::MultiIndex > &globalIndices)
Get the global indices for the provided finite element.
Definition: fehelper.hh:195
Definition: algorithms.hh:17
auto obtainForcesDueToIDBC(const A &assembler)
A helper function to compute the forces due to inhomogeneous Dirichlet BCs.
Definition: functionhelper.hh:170
void forEachLagrangeNodePosition(const LV &localView, F &&f)
A helper function that helps in traversing over the local coordinates of an element and call a user-d...
Definition: traversal.hh:77
decltype(auto) maybeDeref(T &t)
if T is a pointer type, return the dereferenced value, otherwise return the value itself.
Definition: functionhelper.hh:156
def basis(gv, tree)
Definition: basis.py:10
FE class is a base class for all finite elements.
Definition: febase.hh:79
const GridElement & gridElement() const
Get the grid element associated with the local view.
Definition: febase.hh:145
static constexpr int mydim
Dimension of the geometry.
Definition: fetraits.hh:63
Definition: utils/dirichletvalues.hh:38
Type trait to check if a type is a isSharedPtr.
Definition: traits.hh:136
Definition: utils/concepts.hh:641