terms of the local barycentric coordinates of the current tetrahedron, rather than the global coordinate system. Barycentric tracking works on any mesh, irrespective of mesh quality. Particles do not get "lost", and tracking does not require ad-hoc "corrections" or "rescues" to function robustly, because the calculation of particle-face intersections is unambiguous and reproducible, even at small angles of incidence. Each particle position is defined by topology (i.e. the decomposed tet cell it is in) and geometry (i.e. where it is in the cell). No search operations are needed on restart or reconstruct, unlike when particle positions are stored in the global coordinate system. The particle positions file now contains particles' local coordinates and topology, rather than the global coordinates and cell. This change to the output format is not backwards compatible. Existing cases with Lagrangian data will not restart, but they will still run from time zero without any modification. This change was necessary in order to guarantee that the loaded particle is valid, and therefore fundamentally prevent "loss" and "search-failure" type bugs (e.g., 2517, 2442, 2286, 1836, 1461, 1341, 1097). The tracking functions have also been converted to function in terms of displacement, rather than end position. This helps remove floating point error issues, particularly towards the end of a tracking step. Wall bounded streamlines have been removed. The implementation proved incompatible with the new tracking algorithm. ParaView has a surface LIC plugin which provides equivalent, or better, functionality. Additionally, bug report <https://bugs.openfoam.org/view.php?id=2517> is resolved by this change.
300 lines
8.1 KiB
C++
300 lines
8.1 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM 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.
|
|
|
|
OpenFOAM 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 OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
Class
|
|
Foam::streamLineParticle
|
|
|
|
Description
|
|
Particle class that samples fields as it passes through. Used in streamline
|
|
calculation.
|
|
|
|
SourceFiles
|
|
streamLineParticle.C
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#ifndef streamLineParticle_H
|
|
#define streamLineParticle_H
|
|
|
|
#include "particle.H"
|
|
#include "autoPtr.H"
|
|
#include "interpolation.H"
|
|
#include "vectorList.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
|
|
// Forward declaration of friend functions and operators
|
|
|
|
class streamLineParticle;
|
|
|
|
Ostream& operator<<(Ostream&, const streamLineParticle&);
|
|
|
|
/*---------------------------------------------------------------------------*\
|
|
Class streamLineParticle Declaration
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
class streamLineParticle
|
|
:
|
|
public particle
|
|
{
|
|
public:
|
|
|
|
class trackingData
|
|
:
|
|
public particle::TrackingData<Cloud<streamLineParticle>>
|
|
{
|
|
public:
|
|
|
|
// Public data
|
|
|
|
const PtrList<interpolation<scalar>>& vsInterp_;
|
|
|
|
const PtrList<interpolation<vector>>& vvInterp_;
|
|
|
|
const label UIndex_;
|
|
|
|
const bool trackForward_;
|
|
|
|
const label nSubCycle_;
|
|
|
|
const scalar trackLength_;
|
|
|
|
DynamicList<vectorList>& allPositions_;
|
|
|
|
List<DynamicList<scalarList>>& allScalars_;
|
|
|
|
List<DynamicList<vectorList>>& allVectors_;
|
|
|
|
|
|
// Constructors
|
|
|
|
//- Construct from components
|
|
trackingData
|
|
(
|
|
Cloud<streamLineParticle>& cloud,
|
|
const PtrList<interpolation<scalar>>& vsInterp,
|
|
const PtrList<interpolation<vector>>& vvInterp,
|
|
const label UIndex,
|
|
const bool trackForward,
|
|
const label nSubCycle,
|
|
const scalar trackLength,
|
|
DynamicList<List<point>>& allPositions,
|
|
List<DynamicList<scalarList>>& allScalars,
|
|
List<DynamicList<vectorList>>& allVectors
|
|
)
|
|
:
|
|
particle::TrackingData<Cloud<streamLineParticle>>(cloud),
|
|
vsInterp_(vsInterp),
|
|
vvInterp_(vvInterp),
|
|
UIndex_(UIndex),
|
|
trackForward_(trackForward),
|
|
nSubCycle_(nSubCycle),
|
|
trackLength_(trackLength),
|
|
allPositions_(allPositions),
|
|
allScalars_(allScalars),
|
|
allVectors_(allVectors)
|
|
{}
|
|
};
|
|
|
|
|
|
private:
|
|
|
|
// Private data
|
|
|
|
//- Lifetime of particle. Particle dies when reaches 0.
|
|
label lifeTime_;
|
|
|
|
//- Sampled positions
|
|
DynamicList<point> sampledPositions_;
|
|
|
|
//- Sampled scalars
|
|
List<DynamicList<scalar>> sampledScalars_;
|
|
|
|
//- Sampled vectors
|
|
List<DynamicList<vector>> sampledVectors_;
|
|
|
|
|
|
// Private Member Functions
|
|
|
|
//- Interpolate all quantities; return interpolated velocity.
|
|
vector interpolateFields
|
|
(
|
|
const trackingData&,
|
|
const point&,
|
|
const label celli,
|
|
const label facei
|
|
);
|
|
|
|
|
|
public:
|
|
|
|
// Constructors
|
|
|
|
//- Construct from components
|
|
streamLineParticle
|
|
(
|
|
const polyMesh& c,
|
|
const vector& position,
|
|
const label celli,
|
|
const label lifeTime
|
|
);
|
|
|
|
//- Construct from Istream
|
|
streamLineParticle
|
|
(
|
|
const polyMesh& c,
|
|
Istream& is,
|
|
bool readFields = true
|
|
);
|
|
|
|
//- Construct copy
|
|
streamLineParticle(const streamLineParticle& p);
|
|
|
|
//- Construct and return a clone
|
|
autoPtr<particle> clone() const
|
|
{
|
|
return autoPtr<particle>(new streamLineParticle(*this));
|
|
}
|
|
|
|
//- Factory class to read-construct particles used for parallel transfer
|
|
class iNew
|
|
{
|
|
const polyMesh& mesh_;
|
|
|
|
public:
|
|
|
|
iNew(const polyMesh& mesh)
|
|
:
|
|
mesh_(mesh)
|
|
{}
|
|
|
|
autoPtr<streamLineParticle> operator()(Istream& is) const
|
|
{
|
|
return autoPtr<streamLineParticle>
|
|
(
|
|
new streamLineParticle(mesh_, is, true)
|
|
);
|
|
}
|
|
};
|
|
|
|
|
|
// Member Functions
|
|
|
|
// Tracking
|
|
|
|
//- Track all particles to their end point
|
|
bool move(trackingData&, const scalar);
|
|
|
|
//- Overridable function to handle the particle hitting a patch
|
|
// Executed before other patch-hitting functions
|
|
bool hitPatch
|
|
(
|
|
const polyPatch&,
|
|
trackingData& td,
|
|
const label patchi,
|
|
const scalar trackFraction,
|
|
const tetIndices& tetIs
|
|
);
|
|
|
|
//- Overridable function to handle the particle hitting a wedge
|
|
void hitWedgePatch
|
|
(
|
|
const wedgePolyPatch&,
|
|
trackingData& td
|
|
);
|
|
|
|
//- Overridable function to handle the particle hitting a
|
|
// symmetry plane
|
|
void hitSymmetryPlanePatch
|
|
(
|
|
const symmetryPlanePolyPatch&,
|
|
trackingData& td
|
|
);
|
|
|
|
//- Overridable function to handle the particle hitting a
|
|
// symmetry patch
|
|
void hitSymmetryPatch
|
|
(
|
|
const symmetryPolyPatch&,
|
|
trackingData& td
|
|
);
|
|
|
|
//- Overridable function to handle the particle hitting a cyclic
|
|
void hitCyclicPatch
|
|
(
|
|
const cyclicPolyPatch&,
|
|
trackingData& td
|
|
);
|
|
|
|
//- Overridable function to handle the particle hitting a
|
|
//- processorPatch
|
|
void hitProcessorPatch
|
|
(
|
|
const processorPolyPatch&,
|
|
trackingData& td
|
|
);
|
|
|
|
//- Overridable function to handle the particle hitting a wallPatch
|
|
void hitWallPatch
|
|
(
|
|
const wallPolyPatch&,
|
|
trackingData& td,
|
|
const tetIndices&
|
|
);
|
|
|
|
//- Overridable function to handle the particle hitting a polyPatch
|
|
void hitPatch
|
|
(
|
|
const polyPatch&,
|
|
trackingData& td
|
|
);
|
|
|
|
|
|
// I-O
|
|
|
|
//- Read
|
|
static void readFields(Cloud<streamLineParticle>&);
|
|
|
|
//- Write
|
|
static void writeFields(const Cloud<streamLineParticle>&);
|
|
|
|
|
|
// Ostream Operator
|
|
|
|
friend Ostream& operator<<(Ostream&, const streamLineParticle&);
|
|
};
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
} // End namespace Foam
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
#endif
|
|
|
|
// ************************************************************************* //
|