#pragma once #include #include #include #include #include #include #include namespace AMDiS { /** * \addtogroup operators * @{ **/ /// \brief The main implementation of an operator to be used in a \ref LocalAssembler. /** * 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 LocalContext The Element or Intersection type * \tparam GridFunction The GridFunction, a LocalFunction is created from, and * that is evaluated at quadrature points. * * **Requirements:** * - `LocalContext` models either Entity (of codim 0) or Intersection * - `GridFunction` models the \ref Concepts::GridFunction **/ template class GridFunctionOperatorBase : public LocalOperator { using LocalFunction = decltype(localFunction(std::declval())); using Element = typename Impl::ContextTypes::Entity; using Geometry = typename Element::Geometry; using LocalCoordinate = typename GridFunction::EntitySet::LocalCoordinate; using QuadFactory = QuadratureFactory; public: /// \brief Constructor. Stores a copy of `gridFct`. /** * An expression operator takes an expression, following the interface of * \ref ExpressionBase, and stores a copy. Additionally, it gets the * differentiation order, to calculate the quadrature degree in \ref getDegree. **/ GridFunctionOperatorBase(GridFunction const& gridFct, int termOrder) : gridFct_(gridFct) , localFct_(localFunction(gridFct_)) , termOrder_(termOrder) {} // Copy constructor GridFunctionOperatorBase(GridFunctionOperatorBase const& that) : gridFct_(that.gridFct_) , localFct_(localFunction(gridFct_)) // recreate localFct_ on new gridFct_ , termOrder_(that.termOrder_) {} // Move constructor GridFunctionOperatorBase(GridFunctionOperatorBase&& that) : gridFct_(std::move(that.gridFct_)) , localFct_(localFunction(gridFct_)) // recreate localFct_ on new gridFct_ , termOrder_(std::move(that.termOrder_)) {} /// \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`. **/ void bind_impl(Element const& element, Geometry const& geometry) { assert( bool(quadFactory_) ); localFct_.bind(element); quadFactory_->bind(localFct_); } /// Unbinds operator from element. void unbind_impl() { localFct_.unbind(); } template void setQuadFactory(PreQuadFactory const& pre) { quadFactory_ = makeQuadratureFactoryPtr(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 QuadratureCreator 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: /// The gridFunction to be used within the operator GridFunction gridFct_; /// localFunction associated with gridFunction. Mus be updated whenever gridFunction changes. LocalFunction localFct_; /// Assign each element type a quadrature rule std::shared_ptr quadFactory_ = std::make_shared(); /// 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 using Constructable = decltype( new T(std::declval()...) ); public: template ::value, int> = 0> GridFunctionOperatorTransposed(Args&&... args) : transposedOp_(std::forward(args)...) {} template void bind_impl(Element const& element, Geometry const& geometry) { transposedOp_.bind_impl(element, geometry); } void unbind_impl() { transposedOp_.unbind_impl(); } template void setQuadFactory(PreQuadFactory const& pre) { transposedOp_.setQuadFactory(pre); } template void getElementMatrix(Context const& context, RowNode const& rowNode, ColNode const& colNode, ElementMatrix& elementMatrix) { auto elementMatrixTransposed = trans(elementMatrix); transposedOp_.getElementMatrix( context, colNode, rowNode, elementMatrixTransposed); } private: Transposed transposedOp_; }; template struct PreGridFunctionOperator { Tag tag; PreGridFct expr; PreQuadFactory quadFactory; }; /// Store tag and expression to create a \ref GridFunctionOperator template auto makeOperator(Tag tag, PreGridFct const& expr, QuadratureArgs&&... args) { auto preQuadFactory = makePreQuadratureFactory(std::forward(args)...); return PreGridFunctionOperator{tag, expr, preQuadFactory}; } /** @} **/ #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 LocalContext 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, LocalContext, GridFct> {}; /// Generate an \ref GridFunctionOperator from a PreOperator (tag, expr). /// @{ template auto makeLocalOperator(PreGridFunctionOperator const& op, GridView const& gridView) { auto gridFct = makeGridFunction(op.expr, gridView); using GridFctOp = GridFunctionOperator; GridFctOp localOperator{op.tag, gridFct}; localOperator.setQuadFactory(op.quadFactory); return localOperator; } template auto makeLocalOperator(std::reference_wrapper> op, GridView const& gridView) { PreGridFunctionOperator const& op_ref = op; auto gridFct = makeGridFunction(std::ref(op_ref.expr), gridView); using GridFctOp = GridFunctionOperator; GridFctOp localOperator{op_ref.tag, gridFct}; localOperator.setQuadFactory(op_ref.quadFactory); return localOperator; } /// @} /// Generate a shared_ptr to \ref GridFunctionOperator from a PreOperator (tag, expr). /// @{ template auto makeLocalOperatorPtr(PreGridFunctionOperator const& op, GridView const& gridView) { auto gridFct = makeGridFunction(op.expr, gridView); using GridFctOp = GridFunctionOperator; auto localOperator = std::make_shared(op.tag, gridFct); localOperator->setQuadFactory(op.quadFactory); return localOperator; } template auto makeLocalOperatorPtr(std::reference_wrapper> op, GridView const& gridView) { PreGridFunctionOperator const& op_ref = op; auto gridFct = makeGridFunction(std::ref(op_ref.expr), gridView); using GridFctOp = GridFunctionOperator; auto localOperator = std::make_shared(op_ref.tag, gridFct); localOperator->setQuadFactory(op_ref.quadFactory); return localOperator; } /// @} #endif // DOXYGEN } // end namespace AMDiS