Commit a7104b56 authored by Felix Gruber's avatar Felix Gruber
Browse files

Simplify creation of bilinear forms, etc.

The creation of bilinear forms, inner products and linear forms was
always a bit ugly. The interface forced the user to explicitly create
tuples which are mostly an implementation detail of how the terms in
each form are stored.

This commit adds a new factory based interface in which one first has to
specify the underlying space(s). Afterwards one can add terms to the
factory object, one by one. Finally the (bi)linear form or inner product
is generated with a call to the create() method.

This is illustrated on the example of creating a bilinear form.
Previously, this was done as follows:

```
auto bilinearForm = make_BilinearForm(testSpaces, solutionSpaces,
        make_tuple(
            make_IntegralTerm<0,0,IntegrationType::valueValue,
                                  DomainOfIntegration::interior>(cFunc),
            make_IntegralTerm<0,0,IntegrationType::gradValue,
                                  DomainOfIntegration::interior>
                              (minusOneFunc, betaFunc),
            make_IntegralTerm<0,1,IntegrationType::normalVector,
                                  DomainOfIntegration::face>
                              (oneFunc, betaFunc)));
```

With the new interface the same bilinear form can be created with

```
auto bilinearForm
  = bilinearFormWithSpaces(testSpaces, solutionSpaces)
    .addIntegralTerm<0,0,IntegrationType::valueValue,
                         DomainOfIntegration::interior>(c)
    .addIntegralTerm<0,0,IntegrationType::gradValue,
                         DomainOfIntegration::interior>(-1., beta)
    .addIntegralTerm<0,1,IntegrationType::normalVector,
                         DomainOfIntegration::face>(1., beta)
    .create();
```

As one can see there is no need anymore to create
`ConstantGridViewFunction`s for constant coefficients. This is now done
internally by the factory object. Adding terms is done with the
`addIntegralTerm` method which closely follows the `make_IntegralTerm`
function. For LinearTerms, the factory additionally has the methods
`addFunctionalTerm` and `addSkeletalFunctionalTerm`.

Factories have to be initialized with spaces which can easily be done
with the functions
`bilinearFormWithSpaces`,
`innerProductWithSpace` and
`linearFormWithSpace`.

All example programs have been converted to use the new factory classes.
If you want to use them in your own code, remember to include the
headers
`dune/dpg/bilinearformfactory.hh`,
`dune/dpg/innerproductfactory.hh` and
`dune/dpg/linearformfactory.hh`.
parents bf53b8db 9e424b6d
Loading
Loading
Loading
Loading
+3 −0
Original line number Diff line number Diff line
@@ -7,6 +7,7 @@ install(FILES applylocalfunctional_ru_impl.hh
              assemble_helper.hh
              assemble_types.hh
              bilinearform.hh
              bilinearformfactory.hh
              boundarytools.hh
              cholesky.hh
              dpg_system_assembler.hh
@@ -15,12 +16,14 @@ install(FILES applylocalfunctional_ru_impl.hh
              faces.hh
              functionplotter.hh
              innerproduct.hh
              innerproductfactory.hh
              integralterm.hh
              integralterm_rr_impl.hh
              integralterm_ru_impl.hh
              integralterm_uu_impl.hh
              leastsquares.hh
              linearform.hh
              linearformfactory.hh
              linearfunctionalterm.hh
              linearintegralterm.hh
              localcoefficients.hh
+165 −0
Original line number Diff line number Diff line
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_DPG_BILINEARFORMFACTORY_HH
#define DUNE_DPG_BILINEARFORMFACTORY_HH

#include <tuple>
#include <type_traits>

#include <dune/dpg/assemble_types.hh>
#include <dune/dpg/bilinearform.hh>
#include <dune/dpg/functions/gridviewfunctions.hh>
#include <dune/dpg/integralterm.hh>
#include <dune/dpg/spacetuple.hh>

namespace Dune {
  template<class TSpaces, class SolSpaces, class... BilinearTerms>
  class BilinearFormFactory;

  template<class TSpaces, class SolSpaces>
  BilinearFormFactory<TSpaces, SolSpaces>
  bilinearFormWithSpaces(const TSpaces& testSpaces,
                         const SolSpaces& solutionSpaces);

  template<class TSpaces, class SolSpaces, class... BilinearTerms>
  class BilinearFormFactory {
    static_assert(is_SpaceTuplePtr<TSpaces>::value,
        "TSpaces needs to be a SpaceTuplePtr!");
    static_assert(is_SpaceTuplePtr<SolSpaces>::value,
        "SolSpaces needs to be a SpaceTuplePtr!");
    using TestSpacesPtr = TSpaces;
    using SolutionSpacesPtr = SolSpaces;
    using TestSpaces = typename TestSpacesPtr::element_type;
    using SolutionSpaces = typename SolutionSpacesPtr::element_type;
    using BilinearTermsTuple = std::tuple<BilinearTerms...>;
    using GridView = typename std::tuple_element_t<0,TestSpaces>::GridView;

    BilinearFormFactory(const TestSpacesPtr& testSpaces,
                        const SolutionSpacesPtr& solutionSpaces,
                        const BilinearTermsTuple& terms)
      : testSpaces(testSpaces)
      , solutionSpaces(solutionSpaces)
      , terms(terms)
      {}

    friend
    BilinearFormFactory<TSpaces, SolSpaces>
    bilinearFormWithSpaces<TSpaces, SolSpaces>(const TSpaces& testSpaces,
                           const SolSpaces& solutionSpaces);

    template<class TSpaces_, class SolSpaces_, class... BilinearTerms_>
    friend class BilinearFormFactory;

    TestSpacesPtr testSpaces;
    SolutionSpacesPtr solutionSpaces;
    BilinearTermsTuple terms;

  public:
    template<size_t lhsSpaceIndex,
             size_t rhsSpaceIndex,
             IntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor,
             typename std::enable_if<
                         integrationType == IntegrationType::valueValue
                      || integrationType == IntegrationType::normalSign>::type*
                    = nullptr
            >
    auto addIntegralTerm(Factor c)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto newTerm = make_IntegralTerm<lhsSpaceIndex,
                                       rhsSpaceIndex,
                                       integrationType,
                                       domainOfIntegration>(std::move(cFunc));
      using NewTerm = decltype(newTerm);
      using NewBilinearFormFactory
        = BilinearFormFactory<TSpaces, SolSpaces, BilinearTerms..., NewTerm>;
      return NewBilinearFormFactory{
               testSpaces, solutionSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    template<size_t lhsSpaceIndex,
             size_t rhsSpaceIndex,
             IntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor, class Direction,
             typename std::enable_if<
                         integrationType == IntegrationType::gradValue
                      || integrationType == IntegrationType::valueGrad
                      || integrationType == IntegrationType::gradGrad
                      || integrationType == IntegrationType::normalVector
                      || integrationType ==
                                IntegrationType::travelDistanceWeighted
                      >::type*
               = nullptr
            >
    auto addIntegralTerm(Factor c, Direction beta)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto betaFunc = Functions::detail::toGridViewFunction<GridView>(beta);
      auto newTerm = make_IntegralTerm<lhsSpaceIndex,
                                       rhsSpaceIndex,
                                       integrationType,
                                       domainOfIntegration>
                                      (std::move(cFunc), std::move(betaFunc));
      using NewTerm = decltype(newTerm);
      using NewBilinearFormFactory
        = BilinearFormFactory<TSpaces, SolSpaces, BilinearTerms..., NewTerm>;
      return NewBilinearFormFactory{
               testSpaces, solutionSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    template<size_t lhsSpaceIndex,
             size_t rhsSpaceIndex,
             IntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor, class Direction,
             typename std::enable_if<
                         integrationType == IntegrationType::gradGrad>::type*
                    = nullptr
            >
    auto addIntegralTerm(Factor c, Direction lhsBeta, Direction rhsBeta)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto lhsBetaFunc
          = Functions::detail::toGridViewFunction<GridView>(lhsBeta);
      auto rhsBetaFunc
          = Functions::detail::toGridViewFunction<GridView>(rhsBeta);
      auto newTerm = make_IntegralTerm<lhsSpaceIndex,
                                       rhsSpaceIndex,
                                       integrationType,
                                       domainOfIntegration>
                                      (std::move(cFunc),
                                       std::move(lhsBetaFunc),
                                       std::move(rhsBetaFunc));
      using NewTerm = decltype(newTerm);
      using NewBilinearFormFactory
        = BilinearFormFactory<TSpaces, SolSpaces, BilinearTerms..., NewTerm>;
      return NewBilinearFormFactory{
               testSpaces, solutionSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    BilinearForm<TSpaces, SolSpaces, std::tuple<BilinearTerms...>>
    create()
    {
      return {testSpaces, solutionSpaces, terms};
    }
  };

  template<class TSpaces, class SolSpaces>
  BilinearFormFactory<TSpaces, SolSpaces>
  bilinearFormWithSpaces(const TSpaces& testSpaces,
                         const SolSpaces& solutionSpaces)
  {
    return {testSpaces, solutionSpaces, std::tuple<>{}};
  }
}

#endif
+41 −0
Original line number Diff line number Diff line
@@ -6,6 +6,7 @@
#include <type_traits>
#include <utility>

#include <dune/common/typetraits.hh>
#include <dune/common/typeutilities.hh>
#include <dune/functions/gridfunctions/gridviewentityset.hh>

@@ -85,6 +86,46 @@ makeConstantGridViewFunction(Constant&& constant, const GridView&)
  return {std::forward<Constant>(constant)};
}

namespace detail {
  template<class T>
  using has_localFunction_t
      = decltype(localFunction(std::declval<const T&>()));

  template<class GridView, class T,
           typename std::enable_if<Std::is_detected_v<has_localFunction_t, T>>
                        ::type* = nullptr>
  T toGridViewFunction(T&& t) {
    return std::forward<T>(t);
  }

  template<class GridView, class T,
           typename std::enable_if<
                    std::is_arithmetic<std::decay_t<T>>::value>
                              ::type* = nullptr >
  ConstantGridViewFunction<std::decay_t<T>, GridView>
  toGridViewFunction(T&& t) {
    return {std::forward<T>(t)};
  }

  template<class GridView, class T, int dim,
           typename std::enable_if<
                    std::is_arithmetic<T>::value>
                              ::type* = nullptr >
  ConstantGridViewFunction<FieldVector<T,dim>, GridView>
  toGridViewFunction(FieldVector<T,dim>&& t) {
    return {std::move(t)};
  }

  template<class GridView, class T, int dim,
           typename std::enable_if<
                    std::is_arithmetic<T>::value>
                              ::type* = nullptr >
  ConstantGridViewFunction<FieldVector<T,dim>, GridView>
  toGridViewFunction(const FieldVector<T,dim>& t) {
    return {t};
  }
}

template<class Range, class LocalDomain, class LocalContext, class F>
class LocalPiecewiseConstantGridViewFunction
{
+155 −0
Original line number Diff line number Diff line
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_DPG_INNERPRODUCTFACTORY_HH
#define DUNE_DPG_INNERPRODUCTFACTORY_HH

#include <tuple>
#include <type_traits>

#include <dune/dpg/assemble_types.hh>
#include <dune/dpg/functions/gridviewfunctions.hh>
#include <dune/dpg/innerproduct.hh>
#include <dune/dpg/integralterm.hh>
#include <dune/dpg/spacetuple.hh>

namespace Dune {
  template<class TSpaces, class... InnerProductTerms>
  class InnerProductFactory;

  template<class TSpaces>
  InnerProductFactory<TSpaces>
  innerProductWithSpace(const TSpaces& testSpaces);

  template<class TSpaces, class... InnerProductTerms>
  class InnerProductFactory {
    static_assert(is_SpaceTuplePtr<TSpaces>::value,
        "TSpaces needs to be a SpaceTuplePtr!");
    using TestSpacesPtr = TSpaces;
    using TestSpaces = typename TestSpacesPtr::element_type;
    using InnerProductTermsTuple = std::tuple<InnerProductTerms...>;
    using GridView = typename std::tuple_element_t<0,TestSpaces>::GridView;

    InnerProductFactory(const TestSpacesPtr& testSpaces,
                        const InnerProductTermsTuple& terms)
      : testSpaces(testSpaces)
      , terms(terms)
      {}

    friend
    InnerProductFactory<TSpaces>
    innerProductWithSpace<TSpaces>(const TSpaces& testSpaces);

    template<class TSpaces_, class... InnerProductTerms_>
    friend class InnerProductFactory;

    TestSpacesPtr testSpaces;
    InnerProductTermsTuple terms;

  public:
    template<size_t lhsSpaceIndex,
             size_t rhsSpaceIndex,
             IntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor,
             typename std::enable_if<
                         integrationType == IntegrationType::valueValue
                      || integrationType == IntegrationType::normalSign>::type*
                    = nullptr
            >
    auto addIntegralTerm(Factor c)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto newTerm = make_IntegralTerm<lhsSpaceIndex,
                                       rhsSpaceIndex,
                                       integrationType,
                                       domainOfIntegration>(std::move(cFunc));
      using NewTerm = decltype(newTerm);
      using NewInnerProductFactory
        = InnerProductFactory<TSpaces, InnerProductTerms..., NewTerm>;
      return NewInnerProductFactory{
               testSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    template<size_t lhsSpaceIndex,
             size_t rhsSpaceIndex,
             IntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor, class Direction,
             typename std::enable_if<
                         integrationType == IntegrationType::gradValue
                      || integrationType == IntegrationType::valueGrad
                      || integrationType == IntegrationType::gradGrad
                      || integrationType == IntegrationType::normalVector
                      || integrationType ==
                                IntegrationType::travelDistanceWeighted
                      >::type*
               = nullptr
            >
    auto addIntegralTerm(Factor c, Direction beta)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto betaFunc = Functions::detail::toGridViewFunction<GridView>(beta);
      auto newTerm = make_IntegralTerm<lhsSpaceIndex,
                                       rhsSpaceIndex,
                                       integrationType,
                                       domainOfIntegration>
                                      (std::move(cFunc), std::move(betaFunc));
      using NewTerm = decltype(newTerm);
      using NewInnerProductFactory
        = InnerProductFactory<TSpaces, InnerProductTerms..., NewTerm>;
      return NewInnerProductFactory{
               testSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    template<size_t lhsSpaceIndex,
             size_t rhsSpaceIndex,
             IntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor, class Direction,
             typename std::enable_if<
                         integrationType == IntegrationType::gradGrad>::type*
                    = nullptr
            >
    auto addIntegralTerm(Factor c, Direction lhsBeta, Direction rhsBeta)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto lhsBetaFunc
          = Functions::detail::toGridViewFunction<GridView>(lhsBeta);
      auto rhsBetaFunc
          = Functions::detail::toGridViewFunction<GridView>(rhsBeta);
      auto newTerm = make_IntegralTerm<lhsSpaceIndex,
                                       rhsSpaceIndex,
                                       integrationType,
                                       domainOfIntegration>
                                      (std::move(cFunc),
                                       std::move(lhsBetaFunc),
                                       std::move(rhsBetaFunc));
      using NewTerm = decltype(newTerm);
      using NewInnerProductFactory
        = InnerProductFactory<TSpaces, InnerProductTerms..., NewTerm>;
      return NewInnerProductFactory{
               testSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    InnerProduct<TSpaces, std::tuple<InnerProductTerms...>>
    create()
    {
      return {testSpaces, terms};
    }
  };

  template<class TSpaces>
  InnerProductFactory<TSpaces>
  innerProductWithSpace(const TSpaces& testSpaces)
  {
    return {testSpaces, std::tuple<>{}};
  }
}

#endif
+157 −0
Original line number Diff line number Diff line
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_DPG_LINEARFORMFACTORY_HH
#define DUNE_DPG_LINEARFORMFACTORY_HH

#include <tuple>
#include <type_traits>

#include <dune/dpg/assemble_types.hh>
#include <dune/dpg/functions/gridviewfunctions.hh>
#include <dune/dpg/linearform.hh>
#include <dune/dpg/linearfunctionalterm.hh>
#include <dune/dpg/linearintegralterm.hh>
#include <dune/dpg/spacetuple.hh>

namespace Dune {
  template<class TSpaces, class... LinearFormTerms>
  class LinearFormFactory;

  template<class TSpaces>
  LinearFormFactory<TSpaces>
  linearFormWithSpace(const TSpaces& testSpaces);

  template<class TSpaces, class... LinearFormTerms>
  class LinearFormFactory {
    static_assert(is_SpaceTuplePtr<TSpaces>::value,
        "TSpaces needs to be a SpaceTuplePtr!");
    using TestSpacesPtr = TSpaces;
    using TestSpaces = typename TestSpacesPtr::element_type;
    using LinearFormTermsTuple = std::tuple<LinearFormTerms...>;
    using GridView = typename std::tuple_element_t<0,TestSpaces>::GridView;

    LinearFormFactory(const TestSpacesPtr& testSpaces,
                        const LinearFormTermsTuple& terms)
      : testSpaces(testSpaces)
      , terms(terms)
      {}

    friend
    LinearFormFactory<TSpaces>
    linearFormWithSpace<TSpaces>(const TSpaces& testSpaces);

    template<class TSpaces_, class... LinearFormTerms_>
    friend class LinearFormFactory;

    TestSpacesPtr testSpaces;
    LinearFormTermsTuple terms;

  public:
    template<size_t spaceIndex,
             LinearIntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor>
    auto addIntegralTerm(Factor c)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto newTerm = make_LinearIntegralTerm<spaceIndex,
                                             integrationType,
                                             domainOfIntegration>
                                                    (std::move(cFunc));
      using NewTerm = decltype(newTerm);
      using NewLinearFormFactory
        = LinearFormFactory<TSpaces, LinearFormTerms..., NewTerm>;
      return NewLinearFormFactory{
               testSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    template<size_t spaceIndex,
             LinearIntegrationType integrationType,
             DomainOfIntegration domainOfIntegration,
             class Factor,
             class Direction>
    auto addIntegralTerm(Factor c, Direction beta)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto betaFunc = Functions::detail::toGridViewFunction<GridView>(beta);
      auto newTerm = make_LinearIntegralTerm<spaceIndex,
                                             integrationType,
                                             domainOfIntegration>
                              (std::move(cFunc), std::move(betaFunc));
      using NewTerm = decltype(newTerm);
      using NewLinearFormFactory
        = LinearFormFactory<TSpaces, LinearFormTerms..., NewTerm>;
      return NewLinearFormFactory{
               testSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    template<size_t spaceIndex,
             class SolutionSpace,
             class FunctionalVector>
    auto addFunctionalTerm(const FunctionalVector& functionalVector,
                           const SolutionSpace& solutionSpace)
    {
      auto newTerm = make_LinearFunctionalTerm<spaceIndex, SolutionSpace,
                                               FunctionalVector>
                                              (functionalVector, solutionSpace);
      using NewTerm = decltype(newTerm);
      using NewLinearFormFactory
        = LinearFormFactory<TSpaces, LinearFormTerms..., NewTerm>;
      return NewLinearFormFactory{
               testSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    template<size_t spaceIndex,
             IntegrationType integrationType,
             class SolutionSpace,
             class FunctionalVector,
             class Factor, class Direction,
             typename std::enable_if<
                         integrationType == IntegrationType::normalVector
                      || integrationType ==
                                  IntegrationType::travelDistanceWeighted
                      >::type*
               = nullptr
            >
    auto addSkeletalFunctionalTerm(
        const FunctionalVector& functionalVector,
        const SolutionSpace& solutionSpace,
        Factor c, Direction beta)
    {
      auto cFunc = Functions::detail::toGridViewFunction<GridView>(c);
      auto betaFunc = Functions::detail::toGridViewFunction<GridView>(beta);
      auto newTerm = make_SkeletalLinearFunctionalTerm
                      <spaceIndex, integrationType>
                      (functionalVector, solutionSpace,
                       std::move(cFunc), std::move(betaFunc));
      using NewTerm = decltype(newTerm);
      using NewLinearFormFactory
        = LinearFormFactory<TSpaces, LinearFormTerms..., NewTerm>;
      return NewLinearFormFactory{
               testSpaces,
               std::tuple_cat(terms, std::make_tuple(std::move(newTerm)))
             };
    }

    LinearForm<TSpaces, std::tuple<LinearFormTerms...>>
    create()
    {
      return {testSpaces, terms};
    }
  };

  template<class TSpaces>
  LinearFormFactory<TSpaces>
  linearFormWithSpace(const TSpaces& testSpaces)
  {
    return {testSpaces, std::tuple<>{}};
  }
}

#endif
Loading