Medusa  1.1
Coordinate Free Mehless Method implementation
mm::PUApproximant< vec_t > Class Template Reference

#include <PUApproximant_fwd.hpp>

Detailed Description

template<typename vec_t>
class mm::PUApproximant< vec_t >

An efficient partition-of-unity based approximation method for gluing the local approximations together into a smooth field.

To join the many different local approximations together, a weighted average of the relevant solutions at each point is computed. This can be done continuously and efficiently using a compactly supported partition of unity.

There are more reasonable choices for \( r \), such as the distance to the closest computation node, not included in the neighborhood. The larger the radii \( r \) the more prominent the averaging effect will be.

To evaluate approximant at \( N_q \) points, we can construct a spatial search structure that supports radius based nearest-neighbor queries (such as a \( k \)-d tree).

For more details see this paper.

See also
SheppardInterpolant

Definition at line 40 of file PUApproximant_fwd.hpp.

Static Public Member Functions

static vec_t::scalar_t weight (typename vec_t::scalar_t r)
 Computes \( C^2 \) smooth weights \( w \in \left [ 0, 1 \right ] \). More...
 
template<typename scalar_t , typename engine_t >
static eigen_vt< scalar_tevaluate (const DomainDiscretization< vec_t > &domain, const eigen_vt< scalar_t > &values, const Range< vec_t > &query_points, d_scalar_t radius_factor, const engine_t &engine)
 Evaluates the Partion-of-unity approximant at points 'query_points' constructed on given \( (domain, values) \) pairs. More...
 
template<typename scalar_t , typename engine_t >
static eigen_vt< scalar_tevaluate (const DomainDiscretization< vec_t > &domain, const eigen_vt< scalar_t > &values, const KDTree< vec_t > &query_points_tree, d_scalar_t radius_factor, const engine_t &engine)
 Overload of evaluate using a precomputed \( k \)-d tree. More...
 

Private Types

template<typename scalar_t >
using eigen_vt = Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 >
 Vector type. More...
 
using d_scalar_t = typename vec_t::scalar_t
 Scalar type. More...
 

Member Function Documentation

◆ evaluate() [1/2]

template<typename vec_t >
template<typename scalar_t , typename engine_t >
static eigen_vt<scalar_t> mm::PUApproximant< vec_t >::evaluate ( const DomainDiscretization< vec_t > &  domain,
const eigen_vt< scalar_t > &  values,
const KDTree< vec_t > &  query_points_tree,
d_scalar_t  radius_factor,
const engine_t &  engine 
)
static

Overload of evaluate using a precomputed \( k \)-d tree.

◆ evaluate() [2/2]

template<typename vec_t >
template<typename scalar_t , typename engine_t >
static eigen_vt<scalar_t> mm::PUApproximant< vec_t >::evaluate ( const DomainDiscretization< vec_t > &  domain,
const eigen_vt< scalar_t > &  values,
const Range< vec_t > &  query_points,
d_scalar_t  radius_factor,
const engine_t &  engine 
)
static

Evaluates the Partion-of-unity approximant at points 'query_points' constructed on given \( (domain, values) \) pairs.

The order of the approximation is equal to the order of local approximations of \( engine \).

Computational complexity is \( O((N + N_q) \log N_q) \), where \( N \) is the number of points in the domain and \( N_q \) is the number of query points.

Usage example:

// Domain.
BoxShape<Vec1d> b(0, 1);
DomainDiscretization<Vec1d> domain = b.discretizeBoundaryWithStep(0.01);
int support_size = 2;
domain.findSupport(FindClosest(support_size));
// Values.
Eigen::VectorXd values(2);
for (int i = 0; i < values.size(); i++) {
values[i] = static_cast<double>(i);
}
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec1d, ScaleToClosest> approx({3}, Monomials<Vec1d>(2));
// Test points.
Range<Vec1d> test_points(1);
test_points[0] = Vec1d(0.5);
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto value = PUApproximant<Vec1d>::evaluate(domain, values, test_points, effective_radius,
approx); // value = 0.5
Parameters
domainDomain with nodes and stencils.
valuesFunction values in each domain node.
query_pointsEvaluation points where the approximated function should be evaluated.
radius_factorLocal subdomains have radius radius_factor*h, where h is the distance to the farthest support node. Radius factor must be large enough that local subdomains cover all query points.
engineLocal approximation engine. If the engine is interpolatory, so is the PU approximant.
Returns
Approximated function values in query points.
Examples
test/interpolants/PUApproximant_test.cpp.

◆ weight()

template<typename vec_t >
static vec_t::scalar_t mm::PUApproximant< vec_t >::weight ( typename vec_t::scalar_t  r)
inlinestatic

Computes \( C^2 \) smooth weights \( w \in \left [ 0, 1 \right ] \).

There are more reasonable choices for \( r \), such as the distance to the closest computation node, not included in the neighborhood. The larger the radii \( r \) the more prominent the averaging effect will be.

Parameters
rEffective radius.
Returns
vec_t::scalar_t

Definition at line 56 of file PUApproximant_fwd.hpp.


The documentation for this class was generated from the following file:
mm::PUApproximant::evaluate
static eigen_vt< scalar_t > evaluate(const DomainDiscretization< vec_t > &domain, const eigen_vt< scalar_t > &values, const Range< vec_t > &query_points, d_scalar_t radius_factor, const engine_t &engine)
Evaluates the Partion-of-unity approximant at points 'query_points' constructed on given pairs.
mm::Vec1d
Vec< double, 1 > Vec1d
Convenience typedef for 1d vector of doubles.
Definition: Vec_fwd.hpp:33