<html><head></head><body><div style="font-family: Verdana;font-size: 12.0px;"><div>
<div>Hi Marco,</div>

<div> </div>

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

<div>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)... </div>

<div> </div>

<div>I guess, you will need a GridFunctionAdapter around your LocalAnalyticalFunction.</div>

<div> </div>

<div>Maybe the following code part is helpful for you?</div>

<div> </div>

<div>
<div>  // create grid part<br/>
  typedef Dune::Fem::AdaptiveLeafGridPart< GridType > GridPartType;<br/>
  GridPartType gridPart( grid );</div>

<div>  // function space type<br/>
  typedef typename GridPartType::ctype DomainFieldType;<br/>
  static const int dimDomain = GridPartType::dimensionworld;<br/>
  typedef Dune::Fem::FunctionSpace< DomainFieldType, double, dimDomain, dimRange > FunctionSpaceType;</div>

<div>  // create exact solution<br/>
  typedef Dune::Fem::ExactSolution< FunctionSpaceType > ExactSolutionType;<br/>
  ExactSolutionType exactSolution;<br/>
  typedef Dune::Fem::GridFunctionAdapter< ExactSolutionType, GridPartType > GridExactSolutionType;<br/>
  GridExactSolutionType gridExactSolution( "exact solution", exactSolution, gridPart, 5 );</div>

<div>  // create discrete function space<br/>
  typedef Dune::Fem::DiscontinuousGalerkinSpace< FunctionSpaceType, GridPartType, polOrder > DiscreteFunctionSpaceType;<br/>
  DiscreteFunctionSpaceType discreteFunctionSpace( gridPart );</div>

<div>  // create discrete function<br/>
  typedef Dune::Fem::AdaptiveDiscreteFunction< DiscreteFunctionSpaceType > DiscreteFunctionType;<br/>
  DiscreteFunctionType solution( "solution", discreteFunctionSpace );<br/>
  solution.clear();</div>

<div>  // perform the interpolation<br/>
  interpolate( gridExactSolution, solution );</div>

<div> </div>
</div>

<div> </div>

<div> </div>

<div>Best regards</div>

<div>Stefan.</div>

<div> 
<div name="quote" style="margin:10px 5px 5px 10px; padding: 10px 0 10px 10px; border-left:2px solid #C3D9E5; word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space;">
<div style="margin:0 0 10px 0;"><b>Gesendet:</b> Donnerstag, 04. Juni 2015 um 19:41 Uhr<br/>
<b>Von:</b> "Agnese, Marco" <m.agnese13@imperial.ac.uk><br/>
<b>An:</b> "dune-fem@dune-project.org" <dune-fem@dune-project.org><br/>
<b>Betreff:</b> [dune-fem] LocalFunctionAdapter</div>

<div name="quoted-content">Hi,<br/>
I am trying (without much success :)) to interpolate a c++ analytical function over a grid.<br/>
<br/>
So let's suppose that I have an identical function<br/>
RangeType f(const DomainType& x,const double& t,const EntityType& entity)<br/>
{<br/>
return x;<br/>
}<br/>
<br/>
I have created the following class<br/>
<br/>
template<typename DiscreteSpaceImp><br/>
class LocalAnalyticalFunction<br/>
{<br/>
public:<br/>
typedef DiscreteSpaceImp DiscreteFunctionSpaceType;<br/>
typedef typename DiscreteFunctionSpaceType::FunctionSpaceType FunctionSpaceType;<br/>
typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;<br/>
typedef typename DiscreteFunctionSpaceType::EntityType EntityType;<br/>
<br/>
typedef typename FunctionSpaceType::DomainType DomainType;<br/>
typedef typename FunctionSpaceType::RangeType RangeType;<br/>
<br/>
typedef std::function<RangeType(const DomainType&,const double&,const EntityType&)> AnalyticalFunctionType;<br/>
<br/>
// constructor<br/>
LocalAnalyticalFunction(const AnalyticalFunctionType& f,const double& t):<br/>
f_(f),t_(t)<br/>
{}<br/>
<br/>
// evaluate local function<br/>
inline void evaluate(const DomainType& x,RangeType& val)<br/>
{<br/>
val=f_(entity().geometry().global(x),t_,entity());<br/>
}<br/>
<br/>
// initialize to new entity<br/>
inline void init(const EntityType& entity)<br/>
{<br/>
entity_=&entity;<br/>
}<br/>
<br/>
inline const EntityType& entity() const<br/>
{<br/>
return *entity_;<br/>
}<br/>
<br/>
private:<br/>
EntityType const * entity_;<br/>
const AnalyticalFunctionType& f_;<br/>
const double& t_;<br/>
};<br/>
<br/>
Now I try to interpolate f() doing this<br/>
<br/>
typedef LocalAnalyticalFunction<DiscreteSpaceType> LocalAnalyticalFunctionType;<br/>
LocalAnalyticalFunctionType localAnalyticalFunction(f,t);<br/>
<br/>
// create adapter<br/>
typedef LocalFunctionAdapter<LocalAnalyticalFunctionType> AdaptedFunctionType;<br/>
AdaptedFunctionType fAdapted("adapted function",localAnalyticalFunction,df.gridPart());<br/>
<br/>
interpolate(fAdapted,df);<br/>
<br/>
where df is a discrete function defined over a LagrangeSpace.<br/>
<br/>
<br/>
The interpolate() function returns me the following error<br/>
<br/>
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>&}<br/>
<br/>
therefore seems that it calls the evaluate() using a lagrange point as x which has a type different from DomainType.<br/>
<br/>
What I am doing wrong? I don't understand.<br/>
<br/>
Thank you very much for your help,<br/>
cheers,<br/>
Marco.<br/>
<br/>
<br/>
<br/>
_______________________________________________<br/>
dune-fem mailing list<br/>
dune-fem@dune-project.org<br/>
<a href="http://lists.dune-project.org/mailman/listinfo/dune-fem" target="_blank">http://lists.dune-project.org/mailman/listinfo/dune-fem</a></div>
</div>
</div>
</div></div></body></html>