searchablePlate.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-2017 OpenFOAM Foundation
9  Copyright (C) 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 "searchablePlate.H"
31 #include "SortableList.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(searchablePlate, 0);
39  (
40  searchableSurface,
41  searchablePlate,
42  dict
43  );
45  (
46  searchableSurface,
47  searchablePlate,
48  dict,
49  plate
50  );
51 }
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 Foam::direction Foam::searchablePlate::calcNormal(const point& span)
57 {
58  direction normalDir = 3;
59 
60  for (direction dir = 0; dir < vector::nComponents; ++dir)
61  {
62  if (span[dir] < 0)
63  {
64  // Negative entry. Flag and exit.
65  normalDir = 3;
66  break;
67  }
68  else if (span[dir] < VSMALL)
69  {
70  if (normalDir != 3)
71  {
72  // Multiple zero entries. Flag and exit.
73  normalDir = 3;
74  break;
75  }
76 
77  normalDir = dir;
78  }
79  }
80 
81  if (normalDir == 3)
82  {
84  << "Span should have two positive and one zero entry: "
85  << span << nl
86  << exit(FatalError);
87  }
88 
89  return normalDir;
90 }
91 
92 
93 // Returns miss or hit with face (always 0)
94 Foam::pointIndexHit Foam::searchablePlate::findNearest
95 (
96  const point& sample,
97  const scalar nearestDistSqr
98 ) const
99 {
100  // For every component direction can be
101  // left of min, right of max or inbetween.
102  // - outside points: project first one x plane (either min().x()
103  // or max().x()), then onto y plane and finally z. You should be left
104  // with intersection point
105  // - inside point: find nearest side (compare to mid point). Project onto
106  // that.
107 
108  // Project point on plane.
109  pointIndexHit info(true, sample, 0);
110  info.rawPoint()[normalDir_] = origin_[normalDir_];
111 
112  // Clip to edges if outside
113  for (direction dir = 0; dir < vector::nComponents; ++dir)
114  {
115  if (dir != normalDir_)
116  {
117  if (info.rawPoint()[dir] < origin_[dir])
118  {
119  info.rawPoint()[dir] = origin_[dir];
120  }
121  else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
122  {
123  info.rawPoint()[dir] = origin_[dir]+span_[dir];
124  }
125  }
126  }
127 
128  // Check if outside. Optimisation: could do some checks on distance already
129  // on components above
130  if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
131  {
132  info.setMiss();
133  info.setIndex(-1);
134  }
135 
136  return info;
137 }
138 
139 
140 Foam::pointIndexHit Foam::searchablePlate::findLine
141 (
142  const point& start,
143  const point& end
144 ) const
145 {
146  pointIndexHit info
147  (
148  true,
149  Zero,
150  0
151  );
152 
153  const vector dir(end-start);
154 
155  if (mag(dir[normalDir_]) < VSMALL)
156  {
157  info.setMiss();
158  info.setIndex(-1);
159  }
160  else
161  {
162  scalar t = (origin_[normalDir_]-start[normalDir_]) / dir[normalDir_];
163 
164  if (t < 0 || t > 1)
165  {
166  info.setMiss();
167  info.setIndex(-1);
168  }
169  else
170  {
171  info.rawPoint() = start+t*dir;
172  info.rawPoint()[normalDir_] = origin_[normalDir_];
173 
174  // Clip to edges
175  for (direction dir = 0; dir < vector::nComponents; ++dir)
176  {
177  if (dir != normalDir_)
178  {
179  if (info.rawPoint()[dir] < origin_[dir])
180  {
181  info.setMiss();
182  info.setIndex(-1);
183  break;
184  }
185  else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
186  {
187  info.setMiss();
188  info.setIndex(-1);
189  break;
190  }
191  }
192  }
193  }
194  }
195 
196  // Debug
197  if (info.hit())
198  {
199  treeBoundBox bb(origin_, origin_+span_);
200  bb.min()[normalDir_] -= 1e-6;
201  bb.max()[normalDir_] += 1e-6;
202 
203  if (!bb.contains(info.hitPoint()))
204  {
206  << "bb:" << bb << endl
207  << "origin_:" << origin_ << endl
208  << "span_:" << span_ << endl
209  << "normalDir_:" << normalDir_ << endl
210  << "hitPoint:" << info.hitPoint()
211  << abort(FatalError);
212  }
213  }
214 
215  return info;
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
220 
221 Foam::searchablePlate::searchablePlate
222 (
223  const IOobject& io,
224  const point& origin,
225  const vector& span
226 )
227 :
228  searchableSurface(io),
229  origin_(origin),
230  span_(span),
231  normalDir_(calcNormal(span_))
232 {
234  << " origin:" << origin_
235  << " origin+span:" << origin_+span_
236  << " normal:" << vector::componentNames[normalDir_] << nl;
237 
238  bounds() = boundBox(origin_, origin_ + span_);
239 }
240 
241 
242 Foam::searchablePlate::searchablePlate
243 (
244  const IOobject& io,
245  const dictionary& dict
246 )
247 :
248  searchablePlate(io, dict.get<point>("origin"), dict.get<vector>("span"))
249 {}
250 
251 
252 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
253 
255 {
256  if (regions_.empty())
257  {
258  regions_.resize(1);
259  regions_.first() = "region0";
260  }
261  return regions_;
262 }
263 
264 
266 {
267  return tmp<pointField>::New(1, origin_ + 0.5*span_);
268 }
269 
270 
272 (
273  pointField& centres,
274  scalarField& radiusSqr
275 ) const
276 {
277  centres.resize(1);
278  radiusSqr.resize(1);
279 
280  centres[0] = origin_ + 0.5*span_;
281  radiusSqr[0] = Foam::magSqr(0.5*span_);
282 
283  // Add a bit to make sure all points are tested inside
284  radiusSqr += Foam::sqr(SMALL);
285 }
286 
287 
289 {
290  auto tpts = tmp<pointField>::New(4, origin_);
291  auto& pts = tpts.ref();
292 
293  pts[2] += span_;
294 
295  if (span_.x() < span_.y() && span_.x() < span_.z())
296  {
297  pts[1].y() += span_.y();
298  pts[3].z() += span_.z();
299  }
300  else if (span_.y() < span_.z())
301  {
302  pts[1].x() += span_.x();
303  pts[3].z() += span_.z();
304  }
305  else
306  {
307  pts[1].x() += span_.x();
308  pts[3].y() += span_.y();
309  }
310 
311  return tpts;
312 }
313 
314 
316 {
317  return bb.overlaps(bounds());
318 }
319 
320 
321 void Foam::searchablePlate::findNearest
322 (
323  const pointField& samples,
324  const scalarField& nearestDistSqr,
325  List<pointIndexHit>& info
326 ) const
327 {
328  info.setSize(samples.size());
329 
330  forAll(samples, i)
331  {
332  info[i] = findNearest(samples[i], nearestDistSqr[i]);
333  }
334 }
335 
336 
337 void Foam::searchablePlate::findLine
338 (
339  const pointField& start,
340  const pointField& end,
341  List<pointIndexHit>& info
342 ) const
343 {
344  info.setSize(start.size());
345 
346  forAll(start, i)
347  {
348  info[i] = findLine(start[i], end[i]);
349  }
350 }
351 
352 
354 (
355  const pointField& start,
356  const pointField& end,
357  List<pointIndexHit>& info
358 ) const
359 {
360  findLine(start, end, info);
361 }
362 
363 
365 (
366  const pointField& start,
367  const pointField& end,
369 ) const
370 {
371  List<pointIndexHit> nearestInfo;
372  findLine(start, end, nearestInfo);
373 
374  info.setSize(start.size());
375  forAll(info, pointi)
376  {
377  if (nearestInfo[pointi].hit())
378  {
379  info[pointi].setSize(1);
380  info[pointi][0] = nearestInfo[pointi];
381  }
382  else
383  {
384  info[pointi].clear();
385  }
386  }
387 }
388 
389 
391 (
392  const List<pointIndexHit>& info,
393  labelList& region
394 ) const
395 {
396  region.setSize(info.size());
397  region = 0;
398 }
399 
400 
402 (
403  const List<pointIndexHit>& info,
404  vectorField& normal
405 ) const
406 {
407  normal.setSize(info.size());
408  normal = Zero;
409  forAll(normal, i)
410  {
411  normal[i][normalDir_] = 1.0;
412  }
413 }
414 
415 
417 (
418  const pointField& points,
419  List<volumeType>& volType
420 ) const
421 {
423  << "Volume type not supported for plate."
424  << exit(FatalError);
425 }
426 
427 
428 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::componentNames
static const char *const componentNames[]
Definition: VectorSpace.H:114
Foam::searchablePlate::findLineAll
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
Definition: searchablePlate.C:365
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::searchablePlate::overlaps
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
Definition: searchablePlate.C:315
Foam::searchablePlate
Searching on finite plate. Plate has to be aligned with coordinate axes. Plate defined as origin and ...
Definition: searchablePlate.H:93
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::searchablePlate::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: searchablePlate.C:288
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::boundBox::overlaps
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
SortableList.H
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::Field< vector >
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
searchablePlate.H
samples
scalarField samples(nIntervals, Zero)
Foam::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
dict
dictionary dict
Definition: searchingEngine.H:14
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:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Vector< scalar >
Foam::List< word >
Foam::searchablePlate::coordinates
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
Definition: searchablePlate.C:265
Foam::searchablePlate::findLineAny
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
Definition: searchablePlate.C:354
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::searchablePlate::regions
virtual const wordList & regions() const
Names of regions.
Definition: searchablePlate.C:254
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::searchablePlate::getNormal
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
Definition: searchablePlate.C:402
Foam::searchablePlate::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Definition: searchablePlate.C:272
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::searchablePlate::getVolumeType
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
Definition: searchablePlate.C:417
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::searchablePlate::getRegion
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
Definition: searchablePlate.C:391