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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wallBoundedStreamLine.H"
27 #include "fvMesh.H"
29 #include "sampledSet.H"
30 #include "faceSet.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(wallBoundedStreamLine, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
43 (
44  const PackedBoolList& isWallPatch,
45  const point& seedPt,
46  const label cellI
47 ) const
48 {
49  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
50 
51  const cell& cFaces = mesh.cells()[cellI];
52 
53  label minFaceI = -1;
54  label minTetPtI = -1;
55  scalar minDistSqr = sqr(GREAT);
56 
57  forAll(cFaces, cFaceI)
58  {
59  label faceI = cFaces[cFaceI];
60 
61  if (isWallPatch[faceI])
62  {
63  const face& f = mesh.faces()[faceI];
64  const label fp0 = mesh.tetBasePtIs()[faceI];
65  const point& basePoint = mesh.points()[f[fp0]];
66 
67  label fp = f.fcIndex(fp0);
68  for (label i = 2; i < f.size(); i++)
69  {
70  const point& thisPoint = mesh.points()[f[fp]];
71  label nextFp = f.fcIndex(fp);
72  const point& nextPoint = mesh.points()[f[nextFp]];
73 
74  const triPointRef tri(basePoint, thisPoint, nextPoint);
75 
76  scalar d2 = magSqr(tri.centre() - seedPt);
77  if (d2 < minDistSqr)
78  {
79  minDistSqr = d2;
80  minFaceI = faceI;
81  minTetPtI = i-1;
82  }
83  fp = nextFp;
84  }
85  }
86  }
87 
88  // Put particle in tet
89  return tetIndices
90  (
91  cellI,
92  minFaceI,
93  minTetPtI,
94  mesh
95  );
96 }
97 
98 
100 {
101  //const Time& runTime = obr_.time();
102  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
103 
104 
105  // Determine the 'wall' patches
106  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
107  // These are the faces that need to be followed
108 
110  PackedBoolList isWallPatch(mesh.nFaces());
111  forAll(boundaryPatch().addressing(), i)
112  {
113  isWallPatch[boundaryPatch().addressing()[i]] = 1;
114  }
115 
116 
117 
118  // Find nearest wall particle for the seedPoints
119  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
120 
123  (
124  mesh,
125  cloudName_,
126  initialParticles
127  );
128 
129  {
130  // Get the seed points
131  // ~~~~~~~~~~~~~~~~~~~
132 
133  const sampledSet& seedPoints = sampledSetPtr_();
134 
135 
136  forAll(seedPoints, i)
137  {
138  const point& seedPt = seedPoints[i];
139  label cellI = seedPoints.cells()[i];
140 
141  tetIndices ids(findNearestTet(isWallPatch, seedPt, cellI));
142 
143  if (ids.face() != -1 && isWallPatch[ids.face()])
144  {
145  //Pout<< "Seeding particle :" << nl
146  // << " seedPt:" << seedPt << nl
147  // << " face :" << ids.face() << nl
148  // << " at :" << mesh.faceCentres()[ids.face()] << nl
149  // << " cell :" << mesh.cellCentres()[ids.cell()] << nl
150  // << endl;
151 
152  particles.addParticle
153  (
155  (
156  mesh,
157  ids.faceTri(mesh).centre(),
158  ids.cell(),
159  ids.face(), // tetFace
160  ids.tetPt(),
161  -1, // not on a mesh edge
162  -1, // not on a diagonal edge
163  lifeTime_ // lifetime
164  )
165  );
166  }
167  else
168  {
169  Pout<< type() << " : ignoring seed " << seedPt
170  << " since not in wall cell" << endl;
171  }
172  }
173  }
174 
175  label nSeeds = returnReduce(particles.size(), sumOp<label>());
176 
177  if (log_) Info<< type() << " : seeded " << nSeeds << " particles" << endl;
178 
179 
180 
181  // Read or lookup fields
186  label UIndex = -1;
187 
189  (
190  nSeeds,
191  UIndex,
192  vsFlds,
193  vsInterp,
194  vvFlds,
195  vvInterp
196  );
197 
198  // Additional particle info
200  (
201  particles,
202  vsInterp,
203  vvInterp,
204  UIndex, // index of U in vvInterp
205  trackForward_, // track in +u direction?
206  trackLength_, // fixed track length
207  isWallPatch, // which faces are to follow
208 
209  allTracks_,
210  allScalars_,
212  );
213 
214 
215  // Set very large dt. Note: cannot use GREAT since 1/GREAT is SMALL
216  // which is a trigger value for the tracking...
217  const scalar trackTime = Foam::sqrt(GREAT);
218 
219  // Track
220  particles.move(td, trackTime);
221 }
222 
223 
224 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
225 
227 (
228  const word& name,
229  const objectRegistry& obr,
230  const dictionary& dict,
231  const bool loadFromFiles
232 )
233 :
234  streamLineBase(name, obr, dict, loadFromFiles)
235 {
236  // Check if the available mesh is an fvMesh otherise deactivate
237  if (setActive<fvMesh>())
238  {
239  read(dict_);
240  }
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
245 
247 {}
248 
249 
250 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
251 
253 {
254  if (active_)
255  {
257 
258  // Make sure that the mesh is trackable
259  if (debug)
260  {
261  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
262 
263  // 1. positive volume decomposition tets
264  faceSet faces(mesh, "lowQualityTetFaces", mesh.nFaces()/100+1);
265  if
266  (
268  (
269  mesh,
271  true,
272  &faces
273  )
274  )
275  {
276  label nFaces = returnReduce(faces.size(), sumOp<label>());
277 
279  << "Found " << nFaces
280  <<" faces with low quality or negative volume "
281  << "decomposition tets. Writing to faceSet " << faces.name()
282  << endl;
283  }
284 
285  // 2. all edges on a cell having two faces
286  EdgeMap<label> numFacesPerEdge;
287  forAll(mesh.cells(), cellI)
288  {
289  const cell& cFaces = mesh.cells()[cellI];
290 
291  numFacesPerEdge.clear();
292 
293  forAll(cFaces, cFaceI)
294  {
295  label faceI = cFaces[cFaceI];
296  const face& f = mesh.faces()[faceI];
297  forAll(f, fp)
298  {
299  const edge e(f[fp], f.nextLabel(fp));
301  numFacesPerEdge.find(e);
302  if (eFnd != numFacesPerEdge.end())
303  {
304  eFnd()++;
305  }
306  else
307  {
308  numFacesPerEdge.insert(e, 1);
309  }
310  }
311  }
312 
313  forAllConstIter(EdgeMap<label>, numFacesPerEdge, iter)
314  {
315  if (iter() != 2)
316  {
318  << "problem cell:" << cellI
319  << abort(FatalError);
320  }
321  }
322  }
323  }
324  }
325 }
326 
327 
328 // ************************************************************************* //
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:64
Foam::Cloud::addParticle
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:162
Foam::streamLineBase::obr_
const objectRegistry & obr_
Database this class is registered to.
Definition: streamLineBase.H:72
Foam::polyMeshTetDecomposition::minTetQuality
static const scalar minTetQuality
Minimum tetrahedron quality.
Definition: polyMeshTetDecomposition.H:64
Foam::streamLineBase::allTracks_
DynamicList< List< point > > allTracks_
All tracks. Per track the points it passed through.
Definition: streamLineBase.H:135
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::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:98
Foam::streamLineBase::allVectors_
List< DynamicList< vectorList > > allVectors_
Per vectorField, per track, the sampled values.
Definition: streamLineBase.H:141
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
wallBoundedStreamLineParticleCloud.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::sampledSet::cells
const labelList & cells() const
Definition: sampledSet.H:263
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::HashTable::insert
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::streamLineBase::log_
Switch log_
Switch to send output to Info as well as to file.
Definition: streamLineBase.H:78
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::streamLineBase::sampledSetPtr_
autoPtr< sampledSet > sampledSetPtr_
Seed set engine.
Definition: streamLineBase.H:120
Foam::streamLineBase::read
virtual void read(const dictionary &)
Read the field average data.
Definition: streamLineBase.C:560
Foam::tetIndices::tetPt
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:60
Foam::wallBoundedStreamLine::findNearestTet
tetIndices findNearestTet(const PackedBoolList &isWallPatch, const point &seedPt, const label cellI) const
Find wall tet on cell.
Definition: wallBoundedStreamLine.C:43
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::triangle::centre
Point centre() const
Return centre (centroid)
Definition: triangleI.H:102
Foam::boundaryPatch
Like polyPatch but without reference to mesh. patchIdentifier::index is not used. Used in boundaryMes...
Definition: boundaryPatch.H:50
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
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:49
Foam::wallBoundedStreamLine::~wallBoundedStreamLine
virtual ~wallBoundedStreamLine()
Destructor.
Definition: wallBoundedStreamLine.C:246
Foam::wallBoundedStreamLine::wallBoundedStreamLine
wallBoundedStreamLine(const wallBoundedStreamLine &)
Disallow default bitwise copy construct.
Foam::wallBoundedStreamLineParticle::trackingData
Class used to pass tracking data to the trackToEdge function.
Definition: wallBoundedStreamLineParticle.H:63
wallBoundedStreamLine.H
sampledSet.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::IDLList
Intrusive doubly-linked list.
Definition: IDLList.H:47
Foam::tetIndices::faceTri
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the.
Definition: tetIndicesI.H:109
Foam::Info
messageStream Info
Foam::wallBoundedStreamLineParticle
Particle class that samples fields as it passes through. Used in streamline calculation.
Definition: wallBoundedStreamLineParticle.H:55
Foam::Cloud::size
label size() const
Definition: Cloud.H:175
faceSet.H
Foam::tetIndices::face
label face() const
Return the face.
Definition: tetIndicesI.H:36
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::streamLineBase::wallPatch
autoPtr< indirectPrimitivePatch > wallPatch() const
Construct patch out of all wall patch faces.
Definition: streamLineBase.C:47
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::HashTable::find
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Foam::tetIndices::cell
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
Foam::streamLineBase::cloudName_
word cloudName_
Optional specified name of particles.
Definition: streamLineBase.H:102
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:73
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::wallBoundedStreamLine::track
virtual void track()
Do the actual tracking to fill the track data.
Definition: wallBoundedStreamLine.C:99
Foam::EdgeMap< label >
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::sumOp
Definition: ops.H:162
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::HashTable::clear
void clear()
Clear all entries from table.
Definition: HashTable.C:473
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::streamLineBase::lifeTime_
label lifeTime_
Maximum lifetime (= number of cells) of particle.
Definition: streamLineBase.H:93
Foam::streamLineBase
Definition: streamLineBase.H:62
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::streamLineBase::trackLength_
scalar trackLength_
Track length.
Definition: streamLineBase.H:96
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::polyMeshTetDecomposition::checkFaceTets
static bool checkFaceTets(const polyMesh &mesh, scalar tol=minTetQuality, const bool report=false, labelHashSet *setPtr=NULL)
Check face-decomposition tet volume.
Definition: polyMeshTetDecomposition.C:367
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::streamLineBase::allScalars_
List< DynamicList< scalarList > > allScalars_
Per scalarField, per track, the sampled values.
Definition: streamLineBase.H:138
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::streamLineBase::trackForward_
bool trackForward_
Whether to use +u or -u.
Definition: streamLineBase.H:90
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::wallBoundedStreamLine::read
virtual void read(const dictionary &)
Read settings.
Definition: wallBoundedStreamLine.C:252
Foam::Cloud::move
void move(TrackData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:187