All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
lineFace.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
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 "lineFace.H"
27 #include "meshSearch.H"
28 #include "DynamicList.H"
29 #include "polyMesh.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace sampledSets
37 {
38  defineTypeNameAndDebug(lineFace, 0);
39  addToRunTimeSelectionTable(sampledSet, lineFace, word);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
45 
46 void Foam::sampledSets::lineFace::calcSamples
47 (
48  const polyMesh& mesh,
49  const meshSearch& searchEngine,
50  const vector& start,
51  const vector& end,
52  DynamicList<point>& samplingPts,
53  DynamicList<label>& samplingCells,
54  DynamicList<label>& samplingFaces,
55  DynamicList<label>& samplingSegments,
56  DynamicList<scalar>& samplingCurveDist
57 )
58 {
59  // Test for an inward pointing first hit. If this exists, then the start
60  // point is outside, and the sampling should begin from the hit.
61  label bHitI = -1;
62  const List<pointIndexHit> bHits =
63  searchEngine.intersections(start, end);
64  if (bHits.size())
65  {
66  const label bFaceI = bHits[0].index();
67  const vector bNormal = normalised(mesh.faceAreas()[bFaceI]);
68  const scalar bDot = bNormal & (end - start);
69  if (bDot < 0)
70  {
71  bHitI = 0;
72  }
73  }
74 
75  // If an initial inward pointing hit was not found then initialise the
76  // within the cell containing the start point
77  vector startPt = vector::max;
78  label startFaceI = -1, startCellI = -1;
79  if (bHitI == -1)
80  {
81  startPt = start;
82  startFaceI = -1;
83  startCellI = searchEngine.findCell(start);
84  }
85 
86  // If neither hits nor a containing cell for the start point were found then
87  // the line is completely outside the mesh. Return without generating any
88  // sampling points.
89  if (bHits.size() == 0 && startCellI == -1)
90  {
91  return;
92  }
93 
94  // If nothing is set by now then something has gone wrong
95  if (bHitI == -1 && startCellI == -1)
96  {
98  << "The initial location for the line " << start << " to " << end
99  << "could not be found" << exit(FatalError);
100  }
101 
102  // Loop over the hits, starting a new segment at every second hit
103  for
104  (
105  label sampleSegmentI = 0;
106  bHitI < bHits.size();
107  bHitI += 2, sampleSegmentI += 1
108  )
109  {
110  // Set the start point and topology, unless starting within a cell
111  if (bHitI != -1)
112  {
113  startPt = bHits[bHitI].hitPoint();
114  startFaceI = bHits[bHitI].index();
115  startCellI = mesh.faceOwner()[startFaceI];
116  }
117 
118  // Create a particle. If we are starting on a boundary face, track
119  // backwards into it so that the particle has the correct topology.
120  passiveParticle sampleParticle(mesh, startPt, startCellI);
121  if (startFaceI != -1)
122  {
123  sampleParticle.track(start - end, 0);
124  if (!sampleParticle.onBoundaryFace())
125  {
127  << "Failed to associate with the starting boundary face"
128  << exit(FatalError);
129  }
130  }
131 
132  // Track until a boundary is hit, appending the face intersections to
133  // the lists of samples
134  while (true)
135  {
136  const point pt = sampleParticle.position();
137  const scalar dist = mag(pt - start);
138  const bool first =
139  samplingSegments.size() == 0
140  || samplingSegments.last() != sampleSegmentI;
141 
142  if (sampleParticle.onFace())
143  {
144  samplingPts.append(pt);
145  samplingCells.append(sampleParticle.cell());
146  samplingFaces.append(sampleParticle.face());
147  samplingSegments.append(sampleSegmentI);
148  samplingCurveDist.append(dist);
149  }
150 
151  const vector s = (1 - dist/mag(end - start))*(end - start);
152 
153  if
154  (
155  (!first && sampleParticle.onBoundaryFace())
156  || sampleParticle.trackToCell(s, 0) == 0
157  )
158  {
159  break;
160  }
161  }
162  }
163 }
164 
165 
166 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
167 
168 void Foam::sampledSets::lineFace::calcSamples
169 (
170  DynamicList<point>& samplingPts,
171  DynamicList<label>& samplingCells,
172  DynamicList<label>& samplingFaces,
173  DynamicList<label>& samplingSegments,
174  DynamicList<scalar>& samplingCurveDist
175 ) const
176 {
177  calcSamples
178  (
179  mesh(),
180  searchEngine(),
181  start_,
182  end_,
183  samplingPts,
184  samplingCells,
185  samplingFaces,
186  samplingSegments,
187  samplingCurveDist
188  );
189 }
190 
191 
192 void Foam::sampledSets::lineFace::genSamples()
193 {
194  DynamicList<point> samplingPts;
195  DynamicList<label> samplingCells;
196  DynamicList<label> samplingFaces;
197  DynamicList<label> samplingSegments;
198  DynamicList<scalar> samplingCurveDist;
199 
200  calcSamples
201  (
202  samplingPts,
203  samplingCells,
204  samplingFaces,
205  samplingSegments,
206  samplingCurveDist
207  );
208 
209  samplingPts.shrink();
210  samplingCells.shrink();
211  samplingFaces.shrink();
212  samplingSegments.shrink();
213  samplingCurveDist.shrink();
214 
215  setSamples
216  (
217  samplingPts,
218  samplingCells,
219  samplingFaces,
220  samplingSegments,
221  samplingCurveDist
222  );
223 }
224 
225 
226 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
227 
229 (
230  const word& name,
231  const polyMesh& mesh,
232  const meshSearch& searchEngine,
233  const dictionary& dict
234 )
235 :
236  sampledSet(name, mesh, searchEngine, dict),
237  start_(dict.lookup("start")),
238  end_(dict.lookup("end"))
239 {
240  genSamples();
241 
242  if (debug)
243  {
244  write(Info);
245  }
246 }
247 
248 
249 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
250 
252 {}
253 
254 
255 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool onBoundaryFace() const
Is the particle on a boundary face?
Definition: particleI.H:271
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
virtual ~lineFace()
Destructor.
Definition: lineFace.C:251
scalar trackToCell(const vector &displacement, const scalar fraction)
As particle::track, but stops when a new cell is reached.
Definition: particle.C:608
Macros for easy insertion into run-time selection tables.
label cell() const
Return current cell particle is in.
Definition: particleI.H:142
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:61
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:160
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:690
A class for handling words, derived from string.
Definition: word.H:59
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1041
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
bool onFace() const
Is the particle on a face?
Definition: particleI.H:259
virtual bool write()
Sample and write.
Definition: sampledSets.C:241
const vectorField & faceAreas() const
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:801
T & last()
Return the last element of the list.
Definition: UListI.H:128
defineTypeNameAndDebug(arcUniform, 0)
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:589
lineFace(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: lineFace.C:229
Copy of base particle.
vector position() const
Return current particle position.
Definition: particleI.H:283
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
T & first()
Return reference to the first element of the list.
Definition: UPtrListI.H:43