mirror of
https://github.com/aportelli/LatAnalyze.git
synced 2024-11-13 01:35:35 +00:00
Merge pull request #3 from mspraggs/master
TabFunction: Added optional interpolation type argument
This commit is contained in:
commit
3f28430e4a
@ -18,6 +18,7 @@ noinst_PROGRAMS = \
|
|||||||
exCompiledDoubleFunction\
|
exCompiledDoubleFunction\
|
||||||
exDerivative \
|
exDerivative \
|
||||||
exIntegrator \
|
exIntegrator \
|
||||||
|
exInterp \
|
||||||
exMat \
|
exMat \
|
||||||
exMathInterpreter \
|
exMathInterpreter \
|
||||||
exPlot \
|
exPlot \
|
||||||
@ -42,6 +43,10 @@ exFit_CFLAGS = -g -O2
|
|||||||
exFit_LDFLAGS = -L../lib/.libs -lLatAnalyze
|
exFit_LDFLAGS = -L../lib/.libs -lLatAnalyze
|
||||||
endif
|
endif
|
||||||
|
|
||||||
|
exInterp_SOURCES = exInterp.cpp
|
||||||
|
exInterp_CFLAGS = -g -O2
|
||||||
|
exInterp_LDFLAGS = -L../lib/.libs -lLatAnalyze
|
||||||
|
|
||||||
exIntegrator_SOURCES = exIntegrator.cpp
|
exIntegrator_SOURCES = exIntegrator.cpp
|
||||||
exIntegrator_CFLAGS = -g -O2
|
exIntegrator_CFLAGS = -g -O2
|
||||||
exIntegrator_LDFLAGS = -L../lib/.libs -lLatAnalyze
|
exIntegrator_LDFLAGS = -L../lib/.libs -lLatAnalyze
|
||||||
|
38
examples/exInterp.cpp
Normal file
38
examples/exInterp.cpp
Normal file
@ -0,0 +1,38 @@
|
|||||||
|
/*
|
||||||
|
* TabFunction.cpp, part of LatAnalyze 3
|
||||||
|
*
|
||||||
|
* Copyright (C) 2013 - 2015 Antonin Portelli, Matt Spraggs
|
||||||
|
*
|
||||||
|
* LatAnalyze 3 is free software: you can redistribute it and/or modify
|
||||||
|
* it under the terms of the GNU General Public License as published by
|
||||||
|
* the Free Software Foundation, either version 3 of the License, or
|
||||||
|
* (at your option) any later version.
|
||||||
|
*
|
||||||
|
* LatAnalyze 3 is distributed in the hope that it will be useful,
|
||||||
|
* but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||||
|
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||||
|
* GNU General Public License for more details.
|
||||||
|
*
|
||||||
|
* You should have received a copy of the GNU General Public License
|
||||||
|
* along with LatAnalyze 3. If not, see <http://www.gnu.org/licenses/>.
|
||||||
|
*/
|
||||||
|
|
||||||
|
#include <LatAnalyze/TabFunction.hpp>
|
||||||
|
|
||||||
|
|
||||||
|
int main(int argc, char* argv[])
|
||||||
|
{
|
||||||
|
Latan::DVec xs(3);
|
||||||
|
xs << -1.0, 0.0, 1.0;
|
||||||
|
Latan::DVec ys(3);
|
||||||
|
ys << 1.0, 0.0, 1.0;
|
||||||
|
|
||||||
|
auto tab = Latan::TabFunction(xs, ys, Latan::InterpType::QUADRATIC);
|
||||||
|
|
||||||
|
std::cout << "Interpolating naive y = x^2 data..." << std::endl;
|
||||||
|
for (double x = -1.0; x < 1.0; x += 0.1) {
|
||||||
|
double y = tab(&x);
|
||||||
|
std::cout << "y @ " << x << " = " << y;
|
||||||
|
std::cout << " ( " << x * x << " expected)" << std::endl;
|
||||||
|
}
|
||||||
|
}
|
@ -27,14 +27,16 @@ using namespace Latan;
|
|||||||
* TabFunction implementation *
|
* TabFunction implementation *
|
||||||
******************************************************************************/
|
******************************************************************************/
|
||||||
// constructors ////////////////////////////////////////////////////////////////
|
// constructors ////////////////////////////////////////////////////////////////
|
||||||
TabFunction::TabFunction(const DVec &x, const DVec &y)
|
TabFunction::TabFunction(const DVec &x, const DVec &y,
|
||||||
: TabFunction()
|
const InterpType interpType)
|
||||||
|
: interpType_(interpType)
|
||||||
{
|
{
|
||||||
setData(x, y);
|
setData(x, y);
|
||||||
}
|
}
|
||||||
|
|
||||||
TabFunction::TabFunction(const XYStatData &data, const Index i, const Index j)
|
TabFunction::TabFunction(const XYStatData &data, const Index i, const Index j,
|
||||||
: TabFunction()
|
const InterpType interpType)
|
||||||
|
: interpType_(interpType)
|
||||||
{
|
{
|
||||||
setData(data, i, j);
|
setData(data, i, j);
|
||||||
}
|
}
|
||||||
@ -60,26 +62,72 @@ void TabFunction::setData(const XYStatData &data, const Index i, const Index j)
|
|||||||
// function call ///////////////////////////////////////////////////////////////
|
// function call ///////////////////////////////////////////////////////////////
|
||||||
double TabFunction::operator()(const double *arg) const
|
double TabFunction::operator()(const double *arg) const
|
||||||
{
|
{
|
||||||
double x = arg[0], x_a, x_b, y_a, y_b;
|
double result = 0.0, x = arg[0];
|
||||||
|
|
||||||
if ((x < value_.begin()->first)||(x >= value_.rbegin()->first))
|
|
||||||
{
|
if ((x < value_.begin()->first) || (x >= value_.rbegin()->first)) {
|
||||||
LATAN_ERROR(Range, "tabulated function variable out of range (x= "
|
LATAN_ERROR(Range, "tabulated function variable out of range "
|
||||||
+ strFrom(x) + " not in ["
|
"(x= " + strFrom(x) + " not in ["
|
||||||
+ strFrom(value_.begin()->first) + ", "
|
+ strFrom(value_.begin()->first) + ", "
|
||||||
+ strFrom(value_.rbegin()->first) + "])");
|
+ strFrom(value_.rbegin()->first) + "])");
|
||||||
|
}
|
||||||
|
|
||||||
|
auto i = value_.equal_range(x);
|
||||||
|
auto low = (x == i.first->first) ? i.first : prev(i.first);
|
||||||
|
auto high = i.second;
|
||||||
|
|
||||||
|
switch (interpType_) {
|
||||||
|
case InterpType::LINEAR: {
|
||||||
|
double x_a, x_b, y_a, y_b;
|
||||||
|
|
||||||
|
x_a = low->first;
|
||||||
|
x_b = high->first;
|
||||||
|
y_a = low->second;
|
||||||
|
y_b = high->second;
|
||||||
|
result = y_a + (x - x_a) * (y_b - y_a) / (x_b - x_a);
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
case InterpType::NEAREST: {
|
||||||
|
result = nearest(x)->second;
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
case InterpType::QUADRATIC: {
|
||||||
|
double xs[3], ys[3], as[3];
|
||||||
|
auto it = nearest(x);
|
||||||
|
if (it == value_.begin()) {
|
||||||
|
it = next(it);
|
||||||
|
}
|
||||||
|
else if (it == prev(value_.end())) {
|
||||||
|
it = prev(it);
|
||||||
|
}
|
||||||
|
xs[0] = prev(it)->first;
|
||||||
|
ys[0] = prev(it)->second;
|
||||||
|
xs[1] = it->first;
|
||||||
|
ys[1] = it->second;
|
||||||
|
xs[2] = next(it)->first;
|
||||||
|
ys[2] = next(it)->second;
|
||||||
|
|
||||||
|
// Lagrange polynomial coefficient computation
|
||||||
|
as[0]
|
||||||
|
= (x - xs[1]) / (xs[0] - xs[1])
|
||||||
|
* (x - xs[2]) / (xs[0] - xs[2]);
|
||||||
|
as[1]
|
||||||
|
= (x - xs[0]) / (xs[1] - xs[0])
|
||||||
|
* (x - xs[2]) / (xs[1] - xs[2]);
|
||||||
|
as[2]
|
||||||
|
= (x - xs[0]) / (xs[2] - xs[0])
|
||||||
|
* (x - xs[1]) / (xs[2] - xs[1]);
|
||||||
|
result = as[0] * ys[0] + as[1] * ys[1] + as[2] * ys[2];
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
default:
|
||||||
|
int intType = static_cast<int>(interpType_);
|
||||||
|
LATAN_ERROR(Implementation, "unsupported interpolation type in "
|
||||||
|
"tabulated function: "
|
||||||
|
+ strFrom(intType));
|
||||||
}
|
}
|
||||||
|
|
||||||
auto i = value_.equal_range(x);
|
return result;
|
||||||
auto low = (x == i.first->first) ? i.first : prev(i.first);
|
|
||||||
auto high = i.second;
|
|
||||||
|
|
||||||
x_a = low->first;
|
|
||||||
x_b = high->first;
|
|
||||||
y_a = low->second;
|
|
||||||
y_b = high->second;
|
|
||||||
|
|
||||||
return y_a + (x - x_a)*(y_b - y_a)/(x_b - x_a);
|
|
||||||
}
|
}
|
||||||
|
|
||||||
// DoubleFunction factory //////////////////////////////////////////////////////
|
// DoubleFunction factory //////////////////////////////////////////////////////
|
||||||
@ -101,13 +149,30 @@ DoubleFunction TabFunction::makeFunction(const bool makeHardCopy) const
|
|||||||
return res;
|
return res;
|
||||||
}
|
}
|
||||||
|
|
||||||
DoubleFunction Latan::interpolate(const DVec &x, const DVec &y)
|
DoubleFunction Latan::interpolate(const DVec &x, const DVec &y,
|
||||||
|
const InterpType interpType)
|
||||||
{
|
{
|
||||||
return TabFunction(x,y).makeFunction();
|
return TabFunction(x, y, interpType).makeFunction();
|
||||||
}
|
}
|
||||||
|
|
||||||
DoubleFunction Latan::interpolate(const XYStatData &data, const Index i,
|
DoubleFunction Latan::interpolate(const XYStatData &data, const Index i,
|
||||||
const Index j)
|
const Index j, const InterpType interpType)
|
||||||
{
|
{
|
||||||
return TabFunction(data, i, j).makeFunction();
|
return TabFunction(data, i, j, interpType).makeFunction();
|
||||||
|
}
|
||||||
|
|
||||||
|
map<double, double>::const_iterator TabFunction::nearest(const double x) const
|
||||||
|
{
|
||||||
|
map<double, double>::const_iterator ret;
|
||||||
|
auto i = value_.equal_range(x);
|
||||||
|
auto low = (x == i.first->first) ? i.first : prev(i.first);
|
||||||
|
auto high = i.second;
|
||||||
|
if (fabs(high->first - x) < fabs(low->first - x)) {
|
||||||
|
ret = high;
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
ret = low;
|
||||||
|
}
|
||||||
|
|
||||||
|
return ret;
|
||||||
}
|
}
|
||||||
|
@ -20,9 +20,11 @@
|
|||||||
#ifndef Latan_TabFunction_hpp_
|
#ifndef Latan_TabFunction_hpp_
|
||||||
#define Latan_TabFunction_hpp_
|
#define Latan_TabFunction_hpp_
|
||||||
|
|
||||||
|
#include <algorithm>
|
||||||
#include <map>
|
#include <map>
|
||||||
#include <LatAnalyze/Global.hpp>
|
#include <LatAnalyze/Global.hpp>
|
||||||
#include <LatAnalyze/Function.hpp>
|
#include <LatAnalyze/Function.hpp>
|
||||||
|
#include <LatAnalyze/Math.hpp>
|
||||||
#include <LatAnalyze/XYStatData.hpp>
|
#include <LatAnalyze/XYStatData.hpp>
|
||||||
|
|
||||||
BEGIN_LATAN_NAMESPACE
|
BEGIN_LATAN_NAMESPACE
|
||||||
@ -31,13 +33,22 @@ BEGIN_LATAN_NAMESPACE
|
|||||||
* tabulated function: 1D only *
|
* tabulated function: 1D only *
|
||||||
******************************************************************************/
|
******************************************************************************/
|
||||||
|
|
||||||
|
enum class InterpType
|
||||||
|
{
|
||||||
|
NEAREST,
|
||||||
|
LINEAR,
|
||||||
|
QUADRATIC
|
||||||
|
};
|
||||||
|
|
||||||
class TabFunction: public DoubleFunctionFactory
|
class TabFunction: public DoubleFunctionFactory
|
||||||
{
|
{
|
||||||
public:
|
public:
|
||||||
// constructors
|
// constructors
|
||||||
TabFunction(void) = default;
|
TabFunction(void) = default;
|
||||||
TabFunction(const DVec &x, const DVec &y);
|
TabFunction(const DVec &x, const DVec &y,
|
||||||
TabFunction(const XYStatData &data, const Index i = 0, const Index j = 0);
|
const InterpType interpType = InterpType::LINEAR);
|
||||||
|
TabFunction(const XYStatData &data, const Index i = 0, const Index j = 0,
|
||||||
|
const InterpType interpType = InterpType::LINEAR);
|
||||||
// destructor
|
// destructor
|
||||||
virtual ~TabFunction(void) = default;
|
virtual ~TabFunction(void) = default;
|
||||||
// access
|
// access
|
||||||
@ -48,12 +59,17 @@ public:
|
|||||||
// factory
|
// factory
|
||||||
virtual DoubleFunction makeFunction(const bool makeHardCopy = true) const;
|
virtual DoubleFunction makeFunction(const bool makeHardCopy = true) const;
|
||||||
private:
|
private:
|
||||||
|
std::map<double, double>::const_iterator nearest(const double x) const;
|
||||||
|
|
||||||
std::map<double, double> value_;
|
std::map<double, double> value_;
|
||||||
|
InterpType interpType_;
|
||||||
};
|
};
|
||||||
|
|
||||||
DoubleFunction interpolate(const DVec &x, const DVec &y);
|
DoubleFunction interpolate(const DVec &x, const DVec &y,
|
||||||
|
const InterpType interpType = InterpType::LINEAR);
|
||||||
DoubleFunction interpolate(const XYStatData &data, const Index i = 0,
|
DoubleFunction interpolate(const XYStatData &data, const Index i = 0,
|
||||||
const Index j = 0);
|
const Index j = 0,
|
||||||
|
const InterpType interpType = InterpType::LINEAR);
|
||||||
|
|
||||||
END_LATAN_NAMESPACE
|
END_LATAN_NAMESPACE
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user