wallBoundedStreamLine.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "wallBoundedStreamLine.H"
31 #include "sampledSet.H"
32 #include "faceSet.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41  defineTypeNameAndDebug(wallBoundedStreamLine, 0);
43  (
44  functionObject,
45  wallBoundedStreamLine,
46  dictionary
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
56 (
57  const bitSet& isWallPatch,
58  const point& seedPt,
59  const label celli
60 ) const
61 {
62  const cell& cFaces = mesh_.cells()[celli];
63 
64  label minFacei = -1;
65  label minTetPti = -1;
66  scalar minDistSqr = sqr(GREAT);
67  point nearestPt(GREAT, GREAT, GREAT);
68 
69  for (label facei : cFaces)
70  {
71  if (isWallPatch[facei])
72  {
73  const face& f = mesh_.faces()[facei];
74  const label fp0 = mesh_.tetBasePtIs()[facei];
75  const point& basePoint = mesh_.points()[f[fp0]];
76 
77  label fp = f.fcIndex(fp0);
78  for (label i = 2; i < f.size(); i++)
79  {
80  const point& thisPoint = mesh_.points()[f[fp]];
81  const label nextFp = f.fcIndex(fp);
82  const point& nextPoint = mesh_.points()[f[nextFp]];
83 
84  const triPointRef tri(basePoint, thisPoint, nextPoint);
85 
86  //const scalar d2 = magSqr(tri.centre() - seedPt);
87  const pointHit nearInfo(tri.nearestPoint(seedPt));
88  const scalar d2 = nearInfo.distance();
89  if (d2 < minDistSqr)
90  {
91  nearestPt = nearInfo.rawPoint();
92  minDistSqr = d2;
93  minFacei = facei;
94  minTetPti = i-1;
95  }
96  fp = nextFp;
97  }
98  }
99  }
100 
101  // Return tet and nearest point on wall triangle
103  (
104  tetIndices
105  (
106  celli,
107  minFacei,
108  minTetPti
109  ),
110  nearestPt
111  );
112 }
113 
114 
116 (
117  const triPointRef& tri,
118  const point& pt
119 ) const
120 {
121  //pointHit nearPt(tri.nearestPoint(pt));
122  //if (nearPt.distance() > ROOTSMALL)
123  //{
124  // FatalErrorInFunction << "tri:" << tri
125  // << " seed:" << pt << exit(FatalError);
126  //}
127 
128  return (1.0 - ROOTSMALL)*pt + ROOTSMALL*tri.centre();
129 }
130 
131 
133 {
134  // Determine the 'wall' patches
135  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
136  // These are the faces that need to be followed
137 
139  bitSet isWallPatch(mesh_.nFaces(), boundaryPatch().addressing());
140 
141 
142  // Find nearest wall particle for the seedPoints
143  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
144 
147  (
148  mesh_,
149  cloudName_,
150  initialParticles
151  );
152 
153  {
154  // Get the seed points
155  // ~~~~~~~~~~~~~~~~~~~
156 
157  const sampledSet& seedPoints = sampledSetPoints();
158 
159  forAll(seedPoints, seedi)
160  {
161  const label celli = seedPoints.cells()[seedi];
162 
163  if (celli != -1)
164  {
165  const point& seedPt = seedPoints[seedi];
166  Tuple2<tetIndices, point> nearestId
167  (
168  findNearestTet(isWallPatch, seedPt, celli)
169  );
170  const tetIndices& ids = nearestId.first();
171 
172  if (ids.face() != -1 && isWallPatch[ids.face()])
173  {
174  //Pout<< "Seeding particle :" << nl
175  // << " seedPt:" << seedPt << nl
176  // << " face :" << ids.face() << nl
177  // << " at :" << mesh_.faceCentres()[ids.face()]
178  // << nl
179  // << " cell :" << mesh_.cellCentres()[ids.cell()]
180  // << nl << endl;
181 
182  particles.addParticle
183  (
185  (
186  mesh_,
187  // Perturb seed point to be inside triangle
188  pushIn(ids.faceTri(mesh_), nearestId.second()),
189  ids.cell(),
190  ids.face(), // tetFace
191  ids.tetPt(),
192  -1, // not on a mesh edge
193  -1, // not on a diagonal edge
194  (trackDir_ == trackDirType::FORWARD), // forward?
195  lifeTime_ // lifetime
196  )
197  );
198 
199  if (trackDir_ == trackDirType::BIDIRECTIONAL)
200  {
201  // Additional particle for other half of track
202  particles.addParticle
203  (
205  (
206  mesh_,
207  // Perturb seed point to be inside triangle
208  pushIn(ids.faceTri(mesh_), nearestId.second()),
209  ids.cell(),
210  ids.face(), // tetFace
211  ids.tetPt(),
212  -1, // not on a mesh edge
213  -1, // not on a diagonal edge
214  true, // forward
215  lifeTime_ // lifetime
216  )
217  );
218  }
219  }
220  else
221  {
222  Pout<< type() << " : ignoring seed " << seedPt
223  << " since not in wall cell." << endl;
224  }
225  }
226  }
227  }
228 
229  const label nSeeds = returnReduce(particles.size(), sumOp<label>());
230 
231  Log << type() << " : seeded " << nSeeds << " particles." << endl;
232 
233 
234 
235  // Read or lookup fields
240 
241  label UIndex = -1;
242 
244  (
245  nSeeds,
246  UIndex,
247  vsFlds,
248  vsInterp,
249  vvFlds,
250  vvInterp
251  );
252 
253  // Additional particle info
255  (
256  particles,
257  vsInterp,
258  vvInterp,
259  UIndex, // index of U in vvInterp
260  trackLength_, // fixed track length
261  isWallPatch, // which faces are to follow
262 
263  allTracks_,
264  allScalars_,
266  );
267 
268 
269  // Set very large dt. Note: cannot use GREAT since 1/GREAT is SMALL
270  // which is a trigger value for the tracking...
271  const scalar trackTime = Foam::sqrt(GREAT);
272 
273  // Track
274  particles.move(particles, td, trackTime);
275 }
276 
277 
278 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
279 
281 (
282  const word& name,
283  const Time& runTime,
284  const dictionary& dict
285 )
286 :
288 {
289  read(dict_);
290 }
291 
292 
294 (
295  const word& name,
296  const Time& runTime,
297  const dictionary& dict,
298  const wordList& fieldNames
299 )
300 :
302 {
303  read(dict_);
304 }
305 
306 
307 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
308 
310 {
312  {
313  // Make sure that the mesh is trackable
314  if (debug)
315  {
316  // 1. Positive volume decomposition tets
317  faceSet faces(mesh_, "lowQualityTetFaces", mesh_.nFaces()/100+1);
318 
319  if
320  (
322  (
323  mesh_,
325  true,
326  &faces
327  )
328  )
329  {
330  label nFaces = returnReduce(faces.size(), sumOp<label>());
331 
333  << "Found " << nFaces
334  <<" faces with low quality or negative volume "
335  << "decomposition tets. Writing to faceSet " << faces.name()
336  << endl;
337  }
338 
339  // 2. All edges on a cell having two faces
340  EdgeMap<label> numFacesPerEdge;
341  for (const cell& cFaces : mesh_.cells())
342  {
343  numFacesPerEdge.clear();
344 
345  for (const label facei : cFaces)
346  {
347  const face& f = mesh_.faces()[facei];
348  forAll(f, fp)
349  {
350  const edge e(f[fp], f.nextLabel(fp));
351 
352  ++(numFacesPerEdge(e, 0));
353  }
354  }
355 
356  forAllConstIters(numFacesPerEdge, iter)
357  {
358  if (iter() != 2)
359  {
361  << "problem cell:" << cFaces
362  << abort(FatalError);
363  }
364  }
365  }
366  }
367  }
368 
369  return true;
370 }
371 
372 
373 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:83
Foam::polyMeshTetDecomposition::checkFaceTets
static bool checkFaceTets(const polyMesh &mesh, scalar tol=minTetQuality, const bool report=false, labelHashSet *setPtr=nullptr)
Check face-decomposition tet volume.
Definition: polyMeshTetDecomposition.C:384
Foam::Cloud::addParticle
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:104
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::functionObjects::streamLineBase::trackLength_
scalar trackLength_
Track length.
Definition: streamLineBase.H:112
Foam::functionObjects::streamLineBase::allScalars_
List< DynamicList< scalarList > > allScalars_
Per scalarField, per track, the sampled values.
Definition: streamLineBase.H:145
Foam::polyMeshTetDecomposition::minTetQuality
static const scalar minTetQuality
Minimum tetrahedron quality.
Definition: polyMeshTetDecomposition.H:66
Log
#define Log
Definition: PDRblock.C:35
Foam::functionObjects::wallBoundedStreamLine::wallBoundedStreamLine
wallBoundedStreamLine(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: wallBoundedStreamLine.C:281
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::tetIndices::cell
label cell() const
Return the cell.
Definition: tetIndicesI.H:31
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
wallBoundedStreamLineParticleCloud.H
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::functionObjects::streamLineBase::sampledSetPoints
const sampledSet & sampledSetPoints() const
Demand driven construction of the sampledSet.
Definition: streamLineBase.C:78
Foam::sampledSet::cells
const labelList & cells() const
Definition: sampledSet.H:296
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::tetIndices::face
label face() const
Return the face.
Definition: tetIndicesI.H:43
Foam::PointHit::distance
scalar distance() const noexcept
Return distance to hit.
Definition: PointHit.H:139
Foam::functionObjects::streamLineBase::wallPatch
autoPtr< indirectPrimitivePatch > wallPatch() const
Construct patch out of all wall patch faces.
Definition: streamLineBase.C:98
Foam::faceSet
A list of face labels.
Definition: faceSet.H:51
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::functionObjects::streamLineBase::read
virtual bool read(const dictionary &)
Read the field average data.
Definition: streamLineBase.C:875
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::triangle::centre
Point centre() const
Return centre (centroid)
Definition: triangleI.H:105
Foam::functionObjects::wallBoundedStreamLine::track
virtual void track()
Do the actual tracking to fill the track data.
Definition: wallBoundedStreamLine.C:132
Foam::boundaryPatch
Like polyPatch but without reference to mesh. Used in boundaryMesh to hold data on patches....
Definition: boundaryPatch.H:54
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
Foam::wallBoundedStreamLineParticleCloud
A Cloud of streamLine particles.
Definition: wallBoundedStreamLineParticleCloud.H:52
Foam::wallBoundedStreamLineParticle::trackingData
Class used to pass tracking data to the trackToEdge function.
Definition: wallBoundedStreamLineParticle.H:74
wallBoundedStreamLine.H
sampledSet.H
Foam::tetIndices::faceTri
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:170
Foam::wallBoundedStreamLineParticle
Particle class that samples fields as it passes through. Used in streamline calculation.
Definition: wallBoundedStreamLineParticle.H:66
faceSet.H
Foam::functionObjects::wallBoundedStreamLine::pushIn
point pushIn(const triPointRef &tri, const point &pt) const
Definition: wallBoundedStreamLine.C:116
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::functionObjects::streamLineBase::initInterpolations
void initInterpolations(const label nSeeds, label &UIndex, PtrList< volScalarField > &vsFlds, PtrList< interpolation< scalar >> &vsInterp, PtrList< volVectorField > &vvFlds, PtrList< interpolation< vector >> &vvInterp)
Initialise fields, interpolators and track storage.
Definition: streamLineBase.C:142
Foam::functionObjects::streamLineBase
Definition: streamLineBase.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::ILList
Template class for intrusive linked lists.
Definition: ILList.H:52
Foam::Cloud::move
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:147
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::functionObjects::streamLineBase::allVectors_
List< DynamicList< vectorList > > allVectors_
Per vectorField, per track, the sampled values.
Definition: streamLineBase.H:148
Foam::tetIndices::tetPt
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:55
Foam::functionObjects::streamLineBase::trackDir_
trackDirType trackDir_
Whether to use +u or -u or both.
Definition: streamLineBase.H:106
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::EdgeMap< label >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::HashTable::clear
void clear()
Clear all entries from table.
Definition: HashTable.C:630
f
labelList f(nPoints)
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::functionObjects::streamLineBase::cloudName_
word cloudName_
Optional specified name of particles.
Definition: streamLineBase.H:118
Foam::functionObject::type
virtual const word & type() const =0
Runtime type information.
Foam::Vector< scalar >
Foam::List< word >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::functionObjects::streamLineBase::lifeTime_
label lifeTime_
Maximum lifetime (= number of cells) of particle.
Definition: streamLineBase.H:109
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::Tuple2::second
const T2 & second() const noexcept
Return second.
Definition: Tuple2.H:130
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::triangle::nearestPoint
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:671
Foam::functionObjects::wallBoundedStreamLine::findNearestTet
Tuple2< tetIndices, point > findNearestTet(const bitSet &isWallPatch, const point &seedPt, const label celli) const
Find wall tet on cell.
Definition: wallBoundedStreamLine.C:56
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:60
Foam::functionObjects::wallBoundedStreamLine::read
virtual bool read(const dictionary &)
Read settings.
Definition: wallBoundedStreamLine.C:309
Foam::functionObjects::streamLineBase::allTracks_
DynamicList< List< point > > allTracks_
All tracks. Per track the points it passed through.
Definition: streamLineBase.H:142
Foam::Tuple2::first
const T1 & first() const noexcept
Return first.
Definition: Tuple2.H:118
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54