Discontinuous Galerkin Library
#include "dg/algorithm.h"
Loading...
Searching...
No Matches
dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 > Struct Template Reference

Evaluate a least squares fit More...

Public Types

using value_type = get_value_type<ContainerType0>
 
using container_type = ContainerType0
 

Public Member Functions

 LeastSquaresExtrapolation ()
 
 LeastSquaresExtrapolation (unsigned max, const ContainerType0 &copyable0, const ContainerType1 &copyable1)
 Set maximum number of vectors and allocate memory.
 
void set_max (unsigned max, const ContainerType0 &copyable0, const ContainerType1 &copyable1)
 Set maximum number of vectors and allocate memory.
 
unsigned get_max () const
 
void extrapolate (double alpha, const ContainerType0 &x, double beta, ContainerType1 &y) const
 extrapolate value at a new unkown value \( y = \alpha f(x) + \beta y \)
 
void extrapolate (const ContainerType0 &x, ContainerType1 &y) const
 extrapolate value at a new unkown value \( y = f(x) \)
 
void update (const ContainerType0 &x_new, const ContainerType1 &y_new)
 insert a new entry / train the machine learning algorithm
 

Detailed Description

template<class ContainerType0, class ContainerType1>
struct dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >

Evaluate a least squares fit

This class gathers pairs of (features, labels) vectors \( (\vec x_i, \vec y_i)\) and then constructs a guess for \( y\) for given unkown \( \vec x\) by constructing the least squares coefficients

\[ \min ||a_i \vec x_i - \vec x||\]

to get

\[ \vec y = a_i \vec y_i\]

Note
In the context of generating initial guesses for a matrix equation from previous solutions this method is equivalent to the "rolling QR" algorithm described in https://arxiv.org/pdf/2009.10863.pdf Austin, A.P. and Chalmers N. and Warburton, T. INITIAL GUESSES FOR SEQUENCES OF LINEAR SYSTEMS IN A GPU-ACCELERATED INCOMPRESSIBLE FLOW SOLVER (2021). This means it works for matrix equations that are constant in time.
This works best if the unkown function \( \vec y = f(\vec x) \) is linear and if the \( x_i\) are orthogonal
//Least-squares fit through many points
dg::Grid1d grid( 0., 2.*M_PI, 1, 200);
// We have a travelling wave at various times
// Note that the xs are linearly dependent e.g. x2 = a*x1 + b*x2
dg::HVec x0 = dg::evaluate( [=](double x) { return sin( x -0.0);}, grid);
dg::HVec x1 = dg::evaluate( [=](double x) { return sin( x -0.1);}, grid);
dg::HVec x2 = dg::evaluate( [=](double x) { return sin( x -0.2);}, grid);
dg::HVec x3 = dg::evaluate( [=](double x) { return sin( x -0.3);}, grid);
// An imaginary associated vector to each xs
dg::HVec y0 = dg::evaluate( [=](double x) { return cos( x -0.0);}, grid);
dg::HVec y1 = dg::evaluate( [=](double x) { return cos( x -0.1);}, grid);
dg::HVec y2 = dg::evaluate( [=](double x) { return cos( x -0.2);}, grid);
dg::HVec y = dg::evaluate( [=](double x) { return cos( x -0.3);}, grid);
dg::HVec guess ( y);
double norm = dg::blas1::dot( y,y);
// Given (x0, y0), (x1, y1), (x2, y2) give a good guess for y at x3
// Since the input xs, ys are linearly dependent we expect the exact solution
extraLS.update( x0, y0);
extraLS.update( x1, y1);
extraLS.extrapolate( x3, guess);
dg::blas1::axpby( 1., y, -1., guess);
double err = sqrt( dg::blas1::dot( guess, guess) /norm);
INFO( "Difference LS 1 is "<<err);
CHECK( err < 1e-13);
extraLS.update( x1, x1); // rejected (already there)
extraLS.update( x2, y2); // also rejected (linearly dependent)
extraLS.extrapolate( x3, guess);
dg::blas1::axpby( 1., y, -1., guess);
err = sqrt( dg::blas1::dot( guess, guess) /norm);
INFO( "Difference LS 2 is "<<err);
CHECK( err < 1e-13);
Template Parameters
ContainerTypeAny class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag. Among others
  • dg::HVec (serial), dg::DVec (cuda / omp), dg::MHVec (mpi + serial) or dg::MDVec (mpi + cuda / omp)
  • std::vector<dg::DVec> (vector of shared device vectors), std::array<double, 4> (array of 4 doubles) or std::map < std::string, dg::DVec> ( a map of named vectors)
  • double (scalar) and other primitive types ...
If there are several ContainerTypes in the argument list, then TensorTraits must exist for all of them
See also
See The dg dispatch system for a detailed explanation of our type dispatch system

Member Typedef Documentation

◆ container_type

template<class ContainerType0 , class ContainerType1 >
using dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::container_type = ContainerType0

◆ value_type

template<class ContainerType0 , class ContainerType1 >
using dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::value_type = get_value_type<ContainerType0>

Constructor & Destructor Documentation

◆ LeastSquaresExtrapolation() [1/2]

template<class ContainerType0 , class ContainerType1 >
dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::LeastSquaresExtrapolation ( )
inline

◆ LeastSquaresExtrapolation() [2/2]

template<class ContainerType0 , class ContainerType1 >
dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::LeastSquaresExtrapolation ( unsigned max,
const ContainerType0 & copyable0,
const ContainerType1 & copyable1 )
inline

Set maximum number of vectors and allocate memory.

Parameters
maxmaximum of vectors to use for fit
copyable0the memory for the x is allocated based on this vector
copyable1the memory for the y is allocated based on this vector

Member Function Documentation

◆ extrapolate() [1/2]

template<class ContainerType0 , class ContainerType1 >
void dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::extrapolate ( const ContainerType0 & x,
ContainerType1 & y ) const
inline

extrapolate value at a new unkown value \( y = f(x) \)

Parameters
x(read only) value to extrapolate for
y(write only) contains extrapolated value on output

◆ extrapolate() [2/2]

template<class ContainerType0 , class ContainerType1 >
void dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::extrapolate ( double alpha,
const ContainerType0 & x,
double beta,
ContainerType1 & y ) const
inline

extrapolate value at a new unkown value \( y = \alpha f(x) + \beta y \)

Parameters
alphaQuality of life parameter
x(read only) value to extrapolate for
betaQuality of life parameter
y(write only) contains extrapolated value on output

◆ get_max()

template<class ContainerType0 , class ContainerType1 >
unsigned dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::get_max ( ) const
inline

return the current extrapolation max This may not coincide with the max set in the constructor if values have not been updated yet

◆ set_max()

template<class ContainerType0 , class ContainerType1 >
void dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::set_max ( unsigned max,
const ContainerType0 & copyable0,
const ContainerType1 & copyable1 )
inline

Set maximum number of vectors and allocate memory.

Parameters
maxmaximum of vectors to use for fit
copyable0the memory for the x is allocated based on this vector
copyable1the memory for the y is allocated based on this vector

◆ update()

template<class ContainerType0 , class ContainerType1 >
void dg::LeastSquaresExtrapolation< ContainerType0, ContainerType1 >::update ( const ContainerType0 & x_new,
const ContainerType1 & y_new )
inline

insert a new entry / train the machine learning algorithm

Parameters
x_newthe input
y_newthe corresponding output
Attention
if x_new is in the span of the existing x (i.e. x_new is a linear combination of the x_i, then the new value pair is rejected and the function returns)
Note
It is very good at recognizing linear dependence because it computes the linear algebra in extended accuracy using exblas

The documentation for this struct was generated from the following file: