#pragma once #include #include #include #include #include #include #include #include namespace AMDiS { /** * \addtogroup operators * @{ **/ /// \brief The main implementation of an CRTP-base class for operators using a grid-function /// coefficient to be used in an \ref Assembler. /** * An Operator that takes a GridFunction as coefficient. * Provides quadrature rules and handles the evaluation of the GridFunction at * local coordinates. * * The class is specialized, by deriving from it, in \ref GridFunctionOperator. * * \tparam Derived The class derived from GridFunctionOperatorBase * \tparam LC The Element or Intersection type * \tparam GF The GridFunction, a LocalFunction is created from, and * that is evaluated at quadrature points. * * **Requirements:** * - `LC` models either Entity (of codim 0) or Intersection * - `GF` models the \ref Concepts::GridFunction **/ template class GridFunctionOperatorBase : public LocalOperator { template friend class LocalOperator; using ContextType = Impl::ContextTypes; using Super = LocalOperator; private: using GridFunction = GF; /// The type of the localFunction associated with the GridFunction using LocalFunction = decltype(localFunction(std::declval())); /// The Codim=0 entity of the grid, the localFunction can be bound to using Element = typename ContextType::Entity; /// The geometry-type of the grid element using Geometry = typename Element::Geometry; /// The type of the local coordinates in the \ref Element using LocalCoordinate = typename GF::EntitySet::LocalCoordinate; /// A factory for QuadratureRules that incooperate the order of the LocalFunction using QuadFactory = QuadratureFactory; public: /// \brief Constructor. Stores a copy of `gridFct`. /** * A GridFunctionOperator takes a gridFunction and the * differentiation order of the operator, to calculate the * quadrature degree in \ref getDegree. **/ template GridFunctionOperatorBase(GridFct&& gridFct, int termOrder) : gridFct_(FWD(gridFct)) , termOrder_(termOrder) {} /// Create a quadrature factory from a PreQuadratureFactory, e.g. class derived from \ref QuadratureFactory /// \tparam PQF A PreQuadratureFactory template void setQuadFactory(PQF&& pre) { using ctype = typename Geometry::ctype; quadFactory_.emplace( makeQuadratureFactory(FWD(pre))); } protected: /// Return expression value at LocalCoordinates auto coefficient(LocalCoordinate const& local) const { assert( this->bound_ ); return (*localFct_)(local); } /// Create a quadrature rule using the \ref QuadratureFactory by calculating the /// quadrature order necessary to integrate the (bi)linear-form accurately. template auto const& getQuadratureRule(Dune::GeometryType type, Nodes const&... nodes) const { assert( bool(quadFactory_) ); int quadDegree = this->getDegree(termOrder_, quadFactory_->order(), nodes...); return quadFactory_->rule(type, quadDegree); } private: /// \brief Binds operator to `element` and `geometry`. /** * Binding an operator to the currently visited element in grid traversal. * Since all operators need geometry information, the `Element::Geometry` * object `geometry` is created once, during grid traversal, and passed in. * * By default, it binds the \ref localFct_ to the `element` and the Quadrature * factory to the localFunction. **/ void bind_impl(Element const& element, Geometry const& geometry) { assert( bool(quadFactory_) ); localFct_.emplace(localFunction(gridFct_)); localFct_->bind(element); quadFactory_->bind(localFct_.value()); } /// Unbinds operator from element. void unbind_impl() { localFct_->unbind(); localFct_.reset(); } private: /// The gridFunction to be used within the operator GridFunction gridFct_; /// localFunction associated with gridFunction. Mus be updated whenever gridFunction changes. Dune::Std::optional localFct_; /// Assign each element type a quadrature rule Dune::Std::optional quadFactory_; /// the derivative order of this operator (in {0, 1, 2}) int termOrder_ = 0; }; /// \brief The transposed operator, implemented in term of its transposed by /// calling \ref getElementMatrix with inverted arguments. template class GridFunctionOperatorTransposed : public LocalOperator { template friend class LocalOperator; template using Constructable = decltype( new T(std::declval()...) ); public: template ::value, int> = 0> GridFunctionOperatorTransposed(Args&&... args) : transposedOp_(FWD(args)...) {} /// Redirects the setQuadFactory call top the transposed operator /// \tparam PQF A PreQuadratureFactory template void setQuadFactory(PQF&& pre) { transposedOp_.setQuadFactory(FWD(pre)); } private: /// Redirects the bind call top the transposed operator template void bind_impl(Element const& element, Geometry const& geometry) { transposedOp_.bind(element, geometry); } /// Redirects the unbind call top the transposed operator void unbind_impl() { transposedOp_.unbind(); } /// Apply the assembling to the transposed elementMatrix with interchanged row-/colNode /** * \tparam CG ContextGeometry * \tparam RN RowNode * \tparam CN ColNode * \tparam Mat ElementMatrix **/ template void getElementMatrix(CG const& contextGeometry, RN const& rowNode, CN const& colNode, Mat& elementMatrix) { auto elementMatrixTransposed = transposed(elementMatrix); transposedOp_.getElementMatrix( contextGeometry, colNode, rowNode, elementMatrixTransposed); } private: Transposed transposedOp_; }; #ifndef DOXYGEN template struct PreGridFunctionOperator { Tag tag; PreGridFct expr; PQF quadFactory; }; #endif /// Store tag and expression into a \ref PreGridFunctionOperator to create a \ref GridFunctionOperator template auto makeOperator(Tag tag, Expr&& expr, QuadratureArgs&&... args) { auto pqf = makePreQuadratureFactory(FWD(args)...); return PreGridFunctionOperator{tag, FWD(expr), std::move(pqf)}; } /** @} **/ #ifndef DOXYGEN /// The base-template for GridFunctionOperators /** * An operator can specialize this class, by deriving from \ref GridFunctionOperatorBase. * With the generic function \ref makeLocalOperator, an instance is created. To * distinguisch different GridFunction operators, a tag can be provided that has no * other effect. * * \tparam Tag An Identifier for this GridFunctionOperator * \tparam LC An Element or Intersection the operator is evaluated on * \tparam GridFct A GridFunction evaluated in local coordinates on the bound element **/ template class GridFunctionOperator : public GridFunctionOperatorBase, LC, GridFct> {}; template auto makeGridFunctionOperator(Tag tag, GF&& gf, QF&& qf) { GridFunctionOperator gfo{tag, FWD(gf)}; gfo.setQuadFactory(FWD(qf)); return gfo; } /// Generate an \ref GridFunctionOperator from a PreOperator (tag, expr). /// @{ template auto makeLocalOperator(PreGridFunctionOperator op, GridView const& gridView) { auto gf = makeGridFunction(std::move(op.expr), gridView); return makeGridFunctionOperator(op.tag, std::move(gf), std::move(op.quadFactory)); } template auto makeLocalOperator(std::reference_wrapper> op_ref, GridView const& gridView) { PreGridFunctionOperator const& op = op_ref; auto gf = makeGridFunction(std::ref(op.expr), gridView); return makeGridFunctionOperator(op.tag, std::move(gf), op.quadFactory); } /// @} #endif // DOXYGEN } // end namespace AMDiS