ENH: define linear interpolations for scalar, complex, vector, tensor...

- vector, tensor versions are defined component-wise
  to avoid intermediates.

  The base version uses the form "(1-t)*a + t*b" without any bounds
  checking (ie, will also extrapolate).
This commit is contained in:
Mark Olesen
2023-01-13 18:15:13 +01:00
parent ba153df8db
commit 0c756cc676
17 changed files with 349 additions and 16 deletions

View File

@ -5,7 +5,7 @@
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2019-2020 OpenCFD Ltd.
Copyright (C) 2019-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -59,9 +59,24 @@ int main(int argc, char *argv[])
<< "complexVector::one : " << complexVector::one << nl
<< nl;
{
const complex a(0, 1);
const complex b(20, 100);
Info<< "lerp of " << a << " : " << b << nl;
for (const double t : { 0.0, 0.5, 1.0, -0.5, 1.5 })
{
Info<< " " << t << " = " << lerp(a, b, t) << nl;
}
}
for (complex c : { complex{1, 0}, complex{1, 2}} )
{
Info<< nl << "complex : " << c << nl;
Info<< nl << "complex : " << c
<< " mag = " << c.magnitude()
<< " norm = " << c.magSqr()
<< nl;
Info<< "sin: " << sin(c) << nl
<< "pow(3.0f): " << pow(c, 3.0f) << nl

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2013 OpenFOAM Foundation
Copyright (C) 2018-2019 OpenCFD Ltd.
Copyright (C) 2018-2023 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -28,8 +28,7 @@ License
#include "dictionary.H"
#include "primitiveEntry.H"
#include "dimensionedScalar.H"
#include "dimensionedTensor.H"
#include "dimensionedTypes.H"
using namespace Foam;
@ -102,6 +101,18 @@ int main(int argc, char *argv[])
}
}
{
dimensionedLabel a("min", dimTime, -10);
dimensionedLabel b("max", dimTime, 100);
scalar t = 0.5;
Info<< "lerp of" << nl
<< " " << a << nl
<< " " << b << nl
<< " = " << lerp(a, b, t) << nl
<< " vs " << ((1-t)*a + t*b) << nl
<< nl;
}
Pout<< "zero scalar (time): " << dimensionedScalar(dimTime) << endl;
Pout<< "zero vector: " << dimensionedVector(dimLength) << endl;

View File

@ -40,6 +40,8 @@ Description
#include "Tuple2.H"
#include "Switch.H"
#include "dictionary.H"
#include "primitiveFields.H"
#include "labelVector.H"
using namespace Foam;
@ -93,6 +95,13 @@ void printValPair(const char* desc, const T1& val1, const T2& val2)
}
// The 'naive' implementation. Not exact at end points
inline scalar naive_lerp(const scalar a, const scalar b, const scalar t)
{
return a + t*(b - a);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class T>
@ -180,6 +189,113 @@ int main(int argc, char *argv[])
{
unsigned nFail = 0;
{
const float a = 1e8f, b = 1.0f;
Info<< "lerp exact: "
<< (a == lerp(a, b, 0.0f)) << " "
<< (b == lerp(a, b, 1.0f)) << nl;
Info<< "naive lerp exact: "
<< (a == naive_lerp(a, b, 0.0f)) << " "
<< (b == naive_lerp(a, b, 1.0f)) << nl;
}
{
const scalar a = 1e24, b = 1.0;
Info<< "lerp exact: "
<< (a == lerp(a, b, 0.0)) << " "
<< (b == lerp(a, b, 1.0)) << nl;
Info<< "naive lerp exact: "
<< (a == naive_lerp(a, b, 0.0)) << " "
<< (b == naive_lerp(a, b, 1.0)) << nl;
}
{
const vector a(vector::uniform(1e24)), b(vector::uniform(1.0));
Info<<"lerp exact: "
<< (a == lerp(a, b, 0.0f)) << " "
<< (b == lerp(a, b, 1.0f)) << nl;
Info<< "lerp: "
<< lerp(vector::uniform(0), vector::uniform(100), 0.5) << nl;
}
{
const lerpOp1<vector> half(0.5);
const vector a(vector::uniform(20));
const vector b(vector::uniform(100));
Info<< "lerp half: "
<< a << " : " << b << " => " << half(a, b) << nl;
}
{
const labelVector a(labelVector::uniform(10000));
const labelVector b(labelVector::uniform(1));
Info<< "lerp (labelVector) = "
<< lerp(a, b, 0.1) << nl;
}
{
const scalar a(0);
const scalar b(100);
Info<< "lerp of " << a << " : " << b << nl;
for (const double t : { 0.0, 0.5, 1.0, -0.5, 1.5 })
{
Info<< " " << t << " = " << lerp(a, b, t) << nl;
}
}
// No yet
#if 0
{
const label a(10000);
const label b(1);
Info<<"lerp (label) = "
<< label(lerp(a, b, 0.1)) << nl;
}
{
const bool a(true);
const bool b(false);
Info<<"lerp (bool) = "
<< (lerp(a, b, 0.5)) << nl;
}
#endif
{
const sphericalTensor a(10), b(20);
Info<<"lerp exact: "
<< (a == lerp(a, b, 0.0f)) << " "
<< (b == lerp(a, b, 1.0f)) << nl;
// Info<< "lerp: "
// << lerp(vector::uniform(0), vector::uniform(100), 0.5) << nl;
}
{
const tensor a(tensor::uniform(1e24));
const tensor b(tensor::uniform(0));
Info<<"lerp exact: "
<< (a == lerp(a, b, 0.0f)) << " "
<< (b == lerp(a, b, 1.0f)) << nl;
// Info<< "lerp: "
// << lerp(vector::uniform(0), vector::uniform(100), 0.5) << nl;
}
{
Info<< nl << "Test Switch parsing:" << nl;
nFail += testParsing