[dune-fem] LocalFunctionAdapter

Agnese, Marco m.agnese13 at imperial.ac.uk
Fri Jun 5 13:00:53 CEST 2015


Hi Stefan,
thank you very much for your help but I don't understand. 

If you look at the LocalFunctionAdapter it says:

" LocalFunctionAdapter wrapped a class with a local evaluate method into a grid function"

therefore I assume that it transforms my LocalAnalyticalFunction into a GridFunction. 
This grid function has the plain evaluate(...) method which simply call the underlay method of my LocalAnalyticalFunction and it has also the evaluateQuadrature(...) which I suppose is the one called by the caching quadrature. 

The GridFunctionAdpater is very similar to LocalFunctionAdapter but it has a global defined evaluate instead of a local one.

The function interpolate takes a GridFunction as first argument it should work exactly the same if I pass to it a LocalFunctionAdapter or a GridFunctionAdapter since both or them are GridFunction.

Therefore I don't understand why I need to wrap around my LocalAnalyticalFunction a GridFunctionAdapeter. 

Thank you very much,
regards,
Marco. 

 
________________________________________
From: Stefan Girke [Stefan.Girke at gmx.de]
Sent: Thursday, June 04, 2015 10:02 PM
To: Agnese, Marco
Cc: dune-fem at dune-project.org
Subject: Aw: [dune-fem] LocalFunctionAdapter

Hi Marco,

"therefore seems that it calls the evaluate() using a lagrange point as x which has a type different from DomainType."
That is true. This is a caching issue. Explanation can be found in quadrature/quadrature.hh:109. There is one evaluate() method for two different kinds of evaluation points: Wrapped quadrature points and "normal" evaluations (DomainType)...

I guess, you will need a GridFunctionAdapter around your LocalAnalyticalFunction.

Maybe the following code part is helpful for you?

  // create grid part
  typedef Dune::Fem::AdaptiveLeafGridPart< GridType > GridPartType;
  GridPartType gridPart( grid );
  // function space type
  typedef typename GridPartType::ctype DomainFieldType;
  static const int dimDomain = GridPartType::dimensionworld;
  typedef Dune::Fem::FunctionSpace< DomainFieldType, double, dimDomain, dimRange > FunctionSpaceType;
  // create exact solution
  typedef Dune::Fem::ExactSolution< FunctionSpaceType > ExactSolutionType;
  ExactSolutionType exactSolution;
  typedef Dune::Fem::GridFunctionAdapter< ExactSolutionType, GridPartType > GridExactSolutionType;
  GridExactSolutionType gridExactSolution( "exact solution", exactSolution, gridPart, 5 );
  // create discrete function space
  typedef Dune::Fem::DiscontinuousGalerkinSpace< FunctionSpaceType, GridPartType, polOrder > DiscreteFunctionSpaceType;
  DiscreteFunctionSpaceType discreteFunctionSpace( gridPart );
  // create discrete function
  typedef Dune::Fem::AdaptiveDiscreteFunction< DiscreteFunctionSpaceType > DiscreteFunctionType;
  DiscreteFunctionType solution( "solution", discreteFunctionSpace );
  solution.clear();
  // perform the interpolation
  interpolate( gridExactSolution, solution );



Best regards
Stefan.

Gesendet: Donnerstag, 04. Juni 2015 um 19:41 Uhr
Von: "Agnese, Marco" <m.agnese13 at imperial.ac.uk>
An: "dune-fem at dune-project.org" <dune-fem at dune-project.org>
Betreff: [dune-fem] LocalFunctionAdapter
Hi,
I am trying (without much success :)) to interpolate a c++ analytical function over a grid.

So let's suppose that I have an identical function
RangeType f(const DomainType& x,const double& t,const EntityType& entity)
{
return x;
}

I have created the following class

template<typename DiscreteSpaceImp>
class LocalAnalyticalFunction
{
public:
typedef DiscreteSpaceImp DiscreteFunctionSpaceType;
typedef typename DiscreteFunctionSpaceType::FunctionSpaceType FunctionSpaceType;
typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
typedef typename DiscreteFunctionSpaceType::EntityType EntityType;

typedef typename FunctionSpaceType::DomainType DomainType;
typedef typename FunctionSpaceType::RangeType RangeType;

typedef std::function<RangeType(const DomainType&,const double&,const EntityType&)> AnalyticalFunctionType;

// constructor
LocalAnalyticalFunction(const AnalyticalFunctionType& f,const double& t):
f_(f),t_(t)
{}

// evaluate local function
inline void evaluate(const DomainType& x,RangeType& val)
{
val=f_(entity().geometry().global(x),t_,entity());
}

// initialize to new entity
inline void init(const EntityType& entity)
{
entity_=&entity;
}

inline const EntityType& entity() const
{
return *entity_;
}

private:
EntityType const * entity_;
const AnalyticalFunctionType& f_;
const double& t_;
};

Now I try to interpolate f() doing this

typedef LocalAnalyticalFunction<DiscreteSpaceType> LocalAnalyticalFunctionType;
LocalAnalyticalFunctionType localAnalyticalFunction(f,t);

// create adapter
typedef LocalFunctionAdapter<LocalAnalyticalFunctionType> AdaptedFunctionType;
AdaptedFunctionType fAdapted("adapted function",localAnalyticalFunction,df.gridPart());

interpolate(fAdapted,df);

where df is a discrete function defined over a LagrangeSpace.


The interpolate() function returns me the following error

no known conversion for argument 1 from ‘const Dune::Fem::QuadraturePointWrapper<Dune::Fem::CachingPointList<Dune::Fem::LeafGridPart<Dune::GeometryGrid<Dune::AlbertaGrid<2, 2>, Dune::Fem::VertexFunction<Dune::AlbertaGrid<2, 2> >, std::allocator<void> > >, 0, Dune::Fem::LagrangePointSetTraits<double, 2, 2u> > >’ to ‘const DomainType& {aka const Dune::FieldVector<double, 2>&}

therefore seems that it calls the evaluate() using a lagrange point as x which has a type different from DomainType.

What I am doing wrong? I don't understand.

Thank you very much for your help,
cheers,
Marco.



_______________________________________________
dune-fem mailing list
dune-fem at dune-project.org
http://lists.dune-project.org/mailman/listinfo/dune-fem



More information about the dune-fem mailing list