mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: streamLine: user selectable interpolation method
This commit is contained in:
@ -0,0 +1,99 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: 2.0.0 |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
FoamFile
|
||||
{
|
||||
version 2.0;
|
||||
format ascii;
|
||||
class dictionary;
|
||||
object controlDict;
|
||||
}
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
application icoFoam;
|
||||
|
||||
startFrom startTime;
|
||||
|
||||
startTime 0;
|
||||
|
||||
stopAt endTime;
|
||||
|
||||
endTime 0.5;
|
||||
|
||||
deltaT 0.005;
|
||||
|
||||
writeControl timeStep;
|
||||
|
||||
writeInterval 20;
|
||||
|
||||
purgeWrite 0;
|
||||
|
||||
writeFormat ascii;
|
||||
|
||||
writePrecision 6;
|
||||
|
||||
writeCompression off;
|
||||
|
||||
timeFormat general;
|
||||
|
||||
timePrecision 6;
|
||||
|
||||
runTimeModifiable true;
|
||||
|
||||
functions
|
||||
{
|
||||
streamLines
|
||||
{
|
||||
type streamLine;
|
||||
|
||||
// Where to load it from (if not already in solver)
|
||||
functionObjectLibs ("libfieldFunctionObjects.so");
|
||||
|
||||
// Output every
|
||||
outputControl outputTime;
|
||||
// outputInterval 10;
|
||||
|
||||
setFormat vtk; //gnuplot; //xmgr; //raw; //jplot;
|
||||
|
||||
// Velocity field to use for tracking.
|
||||
U U;
|
||||
|
||||
// Interpolation method. Default is cellPoint. See sampleDict.
|
||||
//interpolationScheme pointMVC;
|
||||
|
||||
// Tracked forwards (+U) or backwards (-U)
|
||||
trackForward true;
|
||||
|
||||
// Names of fields to sample. Should contain above velocity field!
|
||||
fields (p U);
|
||||
|
||||
// Steps particles can travel before being removed
|
||||
lifeTime 10000;
|
||||
|
||||
// Number of steps per cell (estimate). Set to 1 to disable subcycling.
|
||||
nSubCycle 5;
|
||||
|
||||
// Cloud name to use
|
||||
cloudName particleTracks;
|
||||
|
||||
// Seeding method. See the sampleSets in sampleDict.
|
||||
seedSampleSet uniform; //cloud;//triSurfaceMeshPointSet;
|
||||
|
||||
uniformCoeffs
|
||||
{
|
||||
type uniform;
|
||||
axis x; //distance;
|
||||
|
||||
start (-0.0205 0.0001 0.00001);
|
||||
end (-0.0205 0.0005 0.00001);
|
||||
nPoints 100;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -33,6 +33,7 @@ License
|
||||
#include "sampledSet.H"
|
||||
#include "globalIndex.H"
|
||||
#include "mapDistribute.H"
|
||||
#include "interpolationCellPoint.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
@ -75,9 +76,9 @@ void Foam::streamLine::track()
|
||||
|
||||
// Read or lookup fields
|
||||
PtrList<volScalarField> vsFlds;
|
||||
PtrList<interpolationCellPoint<scalar> > vsInterp;
|
||||
PtrList<interpolation<scalar> > vsInterp;
|
||||
PtrList<volVectorField> vvFlds;
|
||||
PtrList<interpolationCellPoint<vector> > vvInterp;
|
||||
PtrList<interpolation<vector> > vvInterp;
|
||||
|
||||
label UIndex = -1;
|
||||
|
||||
@ -95,13 +96,29 @@ void Foam::streamLine::track()
|
||||
vsInterp.setSize(vsFlds.size());
|
||||
forAll(vsFlds, i)
|
||||
{
|
||||
vsInterp.set(i, new interpolationCellPoint<scalar>(vsFlds[i]));
|
||||
vsInterp.set
|
||||
(
|
||||
i,
|
||||
interpolation<scalar>::New
|
||||
(
|
||||
interpolationScheme_,
|
||||
vsFlds[i]
|
||||
)
|
||||
);
|
||||
}
|
||||
ReadFields(mesh, objects, vvFlds);
|
||||
vvInterp.setSize(vvFlds.size());
|
||||
forAll(vvFlds, i)
|
||||
{
|
||||
vvInterp.set(i, new interpolationCellPoint<vector>(vvFlds[i]));
|
||||
vvInterp.set
|
||||
(
|
||||
i,
|
||||
interpolation<vector>::New
|
||||
(
|
||||
interpolationScheme_,
|
||||
vvFlds[i]
|
||||
)
|
||||
);
|
||||
}
|
||||
}
|
||||
else
|
||||
@ -146,7 +163,12 @@ void Foam::streamLine::track()
|
||||
vsInterp.set
|
||||
(
|
||||
nScalar++,
|
||||
new interpolationCellPoint<scalar>(f)
|
||||
//new interpolationCellPoint<scalar>(f)
|
||||
interpolation<scalar>::New
|
||||
(
|
||||
interpolationScheme_,
|
||||
f
|
||||
)
|
||||
);
|
||||
}
|
||||
else if (mesh.foundObject<volVectorField>(fields_[i]))
|
||||
@ -164,7 +186,12 @@ void Foam::streamLine::track()
|
||||
vvInterp.set
|
||||
(
|
||||
nVector++,
|
||||
new interpolationCellPoint<vector>(f)
|
||||
//new interpolationCellPoint<vector>(f)
|
||||
interpolation<vector>::New
|
||||
(
|
||||
interpolationScheme_,
|
||||
f
|
||||
)
|
||||
);
|
||||
}
|
||||
}
|
||||
@ -307,6 +334,15 @@ void Foam::streamLine::read(const dictionary& dict)
|
||||
nSubCycle_ = 1;
|
||||
}
|
||||
|
||||
interpolationScheme_ = dict.lookupOrDefault
|
||||
(
|
||||
"interpolationScheme",
|
||||
interpolationCellPoint<scalar>::typeName
|
||||
);
|
||||
|
||||
//Info<< typeName << " using interpolation " << interpolationScheme_
|
||||
// << endl;
|
||||
|
||||
cloudName_ = dict.lookupOrDefault<word>("cloudName", "streamLine");
|
||||
dict.lookup("seedSampleSet") >> seedSet_;
|
||||
|
||||
|
||||
@ -86,6 +86,9 @@ class streamLine
|
||||
//- Field to transport particle with
|
||||
word UName_;
|
||||
|
||||
//- Interpolation scheme to use
|
||||
word interpolationScheme_;
|
||||
|
||||
//- Whether to use +u or -u
|
||||
bool trackForward_;
|
||||
|
||||
|
||||
@ -38,7 +38,7 @@ SourceFiles
|
||||
|
||||
#include "particle.H"
|
||||
#include "autoPtr.H"
|
||||
#include "interpolationCellPoint.H"
|
||||
#include "interpolation.H"
|
||||
#include "vectorList.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
@ -68,8 +68,8 @@ public:
|
||||
public:
|
||||
|
||||
|
||||
const PtrList<interpolationCellPoint<scalar> >& vsInterp_;
|
||||
const PtrList<interpolationCellPoint<vector> >& vvInterp_;
|
||||
const PtrList<interpolation<scalar> >& vsInterp_;
|
||||
const PtrList<interpolation<vector> >& vvInterp_;
|
||||
const label UIndex_;
|
||||
const bool trackForward_;
|
||||
const label nSubCycle_;
|
||||
@ -84,8 +84,8 @@ public:
|
||||
trackingData
|
||||
(
|
||||
Cloud<streamLineParticle>& cloud,
|
||||
const PtrList<interpolationCellPoint<scalar> >& vsInterp,
|
||||
const PtrList<interpolationCellPoint<vector> >& vvInterp,
|
||||
const PtrList<interpolation<scalar> >& vsInterp,
|
||||
const PtrList<interpolation<vector> >& vvInterp,
|
||||
const label UIndex,
|
||||
const bool trackForward,
|
||||
const label nSubCycle,
|
||||
|
||||
Reference in New Issue
Block a user