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  {
65  << "Span should have two positive and one zero entry. Now:"
66  << span << exit(FatalError);
67  }
68  else if (span[dir] < VSMALL)
69  {
70  if (normalDir == 3)
71  {
72  normalDir = dir;
73  }
74  else
75  {
76  // Multiple zero entries. Flag and exit.
77  normalDir = 3;
78  break;
79  }
80  }
81  }
82 
83  if (normalDir == 3)
84  {
86  << "Span should have two positive and one zero entry. Now:"
87  << span << exit(FatalError);
88  }
89 
90  return normalDir;
91 }
92 
93 
94 // Returns miss or hit with face (always 0)
95 Foam::pointIndexHit Foam::searchablePlate::findNearest
96 (
97  const point& sample,
98  const scalar nearestDistSqr
99 ) const
100 {
101  // For every component direction can be
102  // left of min, right of max or inbetween.
103  // - outside points: project first one x plane (either min().x()
104  // or max().x()), then onto y plane and finally z. You should be left
105  // with intersection point
106  // - inside point: find nearest side (compare to mid point). Project onto
107  // that.
108 
109  // Project point on plane.
110  pointIndexHit info(true, sample, 0);
111  info.rawPoint()[normalDir_] = origin_[normalDir_];
112 
113  // Clip to edges if outside
114  for (direction dir = 0; dir < vector::nComponents; ++dir)
115  {
116  if (dir != normalDir_)
117  {
118  if (info.rawPoint()[dir] < origin_[dir])
119  {
120  info.rawPoint()[dir] = origin_[dir];
121  }
122  else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
123  {
124  info.rawPoint()[dir] = origin_[dir]+span_[dir];
125  }
126  }
127  }
128 
129  // Check if outside. Optimisation: could do some checks on distance already
130  // on components above
131  if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
132  {
133  info.setMiss();
134  info.setIndex(-1);
135  }
136 
137  return info;
138 }
139 
140 
141 Foam::pointIndexHit Foam::searchablePlate::findLine
142 (
143  const point& start,
144  const point& end
145 ) const
146 {
147  pointIndexHit info
148  (
149  true,
150  Zero,
151  0
152  );
153 
154  const vector dir(end-start);
155 
156  if (mag(dir[normalDir_]) < VSMALL)
157  {
158  info.setMiss();
159  info.setIndex(-1);
160  }
161  else
162  {
163  scalar t = (origin_[normalDir_]-start[normalDir_]) / dir[normalDir_];
164 
165  if (t < 0 || t > 1)
166  {
167  info.setMiss();
168  info.setIndex(-1);
169  }
170  else
171  {
172  info.rawPoint() = start+t*dir;
173  info.rawPoint()[normalDir_] = origin_[normalDir_];
174 
175  // Clip to edges
176  for (direction dir = 0; dir < vector::nComponents; ++dir)
177  {
178  if (dir != normalDir_)
179  {
180  if (info.rawPoint()[dir] < origin_[dir])
181  {
182  info.setMiss();
183  info.setIndex(-1);
184  break;
185  }
186  else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
187  {
188  info.setMiss();
189  info.setIndex(-1);
190  break;
191  }
192  }
193  }
194  }
195  }
196 
197  // Debug
198  if (info.hit())
199  {
200  treeBoundBox bb(origin_, origin_+span_);
201  bb.min()[normalDir_] -= 1e-6;
202  bb.max()[normalDir_] += 1e-6;
203 
204  if (!bb.contains(info.hitPoint()))
205  {
207  << "bb:" << bb << endl
208  << "origin_:" << origin_ << endl
209  << "span_:" << span_ << endl
210  << "normalDir_:" << normalDir_ << endl
211  << "hitPoint:" << info.hitPoint()
212  << abort(FatalError);
213  }
214  }
215 
216  return info;
217 }
218 
219 
220 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
221 
222 Foam::searchablePlate::searchablePlate
223 (
224  const IOobject& io,
225  const point& origin,
226  const vector& span
227 )
228 :
229  searchableSurface(io),
230  origin_(origin),
231  span_(span),
232  normalDir_(calcNormal(span_))
233 {
235  << " origin:" << origin_
236  << " origin+span:" << origin_+span_
237  << " normal:" << vector::componentNames[normalDir_] << nl;
238 
239  bounds() = boundBox(origin_, origin_ + span_);
240 }
241 
242 
243 Foam::searchablePlate::searchablePlate
244 (
245  const IOobject& io,
246  const dictionary& dict
247 )
248 :
249  searchableSurface(io),
250  origin_(dict.get<point>("origin")),
251  span_(dict.get<vector>("span")),
252  normalDir_(calcNormal(span_))
253 {
255  << " origin:" << origin_
256  << " origin+span:" << origin_+span_
257  << " normal:" << vector::componentNames[normalDir_] << nl;
258 
259  bounds() = boundBox(origin_, origin_ + span_);
260 }
261 
262 
263 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
264 
266 {
267  if (regions_.empty())
268  {
269  regions_.setSize(1);
270  regions_[0] = "region0";
271  }
272  return regions_;
273 }
274 
275 
277 {
278  return tmp<pointField>::New(1, origin_ + 0.5*span_);
279 }
280 
281 
283 (
284  pointField& centres,
285  scalarField& radiusSqr
286 ) const
287 {
288  centres.setSize(1);
289  centres[0] = origin_ + 0.5*span_;
290 
291  radiusSqr.setSize(1);
292  radiusSqr[0] = Foam::magSqr(0.5*span_);
293 
294  // Add a bit to make sure all points are tested inside
295  radiusSqr += Foam::sqr(SMALL);
296 }
297 
298 
300 {
301  auto tpts = tmp<pointField>::New(4);
302  auto& pts = tpts.ref();
303 
304  pts[0] = origin_;
305  pts[2] = origin_ + span_;
306 
307  if (span_.x() < span_.y() && span_.x() < span_.z())
308  {
309  pts[1] = origin_ + point(0, span_.y(), 0);
310  pts[3] = origin_ + point(0, 0, span_.z());
311  }
312  else if (span_.y() < span_.z())
313  {
314  pts[1] = origin_ + point(span_.x(), 0, 0);
315  pts[3] = origin_ + point(0, 0, span_.z());
316  }
317  else
318  {
319  pts[1] = origin_ + point(span_.x(), 0, 0);
320  pts[3] = origin_ + point(0, span_.y(), 0);
321  }
322 
323  return tpts;
324 }
325 
326 
328 {
329  return
330  (
331  (origin_.x() + span_.x()) >= bb.min().x()
332  && origin_.x() <= bb.max().x()
333  && (origin_.y() + span_.y()) >= bb.min().y()
334  && origin_.y() <= bb.max().y()
335  && (origin_.z() + span_.z()) >= bb.min().z()
336  && origin_.z() <= bb.max().z()
337  );
338 }
339 
340 
341 void Foam::searchablePlate::findNearest
342 (
343  const pointField& samples,
344  const scalarField& nearestDistSqr,
345  List<pointIndexHit>& info
346 ) const
347 {
348  info.setSize(samples.size());
349 
350  forAll(samples, i)
351  {
352  info[i] = findNearest(samples[i], nearestDistSqr[i]);
353  }
354 }
355 
356 
357 void Foam::searchablePlate::findLine
358 (
359  const pointField& start,
360  const pointField& end,
361  List<pointIndexHit>& info
362 ) const
363 {
364  info.setSize(start.size());
365 
366  forAll(start, i)
367  {
368  info[i] = findLine(start[i], end[i]);
369  }
370 }
371 
372 
374 (
375  const pointField& start,
376  const pointField& end,
377  List<pointIndexHit>& info
378 ) const
379 {
380  findLine(start, end, info);
381 }
382 
383 
385 (
386  const pointField& start,
387  const pointField& end,
389 ) const
390 {
391  List<pointIndexHit> nearestInfo;
392  findLine(start, end, nearestInfo);
393 
394  info.setSize(start.size());
395  forAll(info, pointi)
396  {
397  if (nearestInfo[pointi].hit())
398  {
399  info[pointi].setSize(1);
400  info[pointi][0] = nearestInfo[pointi];
401  }
402  else
403  {
404  info[pointi].clear();
405  }
406  }
407 }
408 
409 
411 (
412  const List<pointIndexHit>& info,
413  labelList& region
414 ) const
415 {
416  region.setSize(info.size());
417  region = 0;
418 }
419 
420 
422 (
423  const List<pointIndexHit>& info,
424  vectorField& normal
425 ) const
426 {
427  normal.setSize(info.size());
428  normal = Zero;
429  forAll(normal, i)
430  {
431  normal[i][normalDir_] = 1.0;
432  }
433 }
434 
435 
437 (
438  const pointField& points,
439  List<volumeType>& volType
440 ) const
441 {
443  << "Volume type not supported for plate."
444  << exit(FatalError);
445 }
446 
447 
448 // ************************************************************************* //
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::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Foam::VectorSpace< Vector< scalar >, scalar, 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:385
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
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::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
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:327
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
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:299
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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:55
Foam::Field< vector >
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:360
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)
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:137
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
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
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:276
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:374
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:47
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:265
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:422
Foam::searchablePlate::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Definition: searchablePlate.C:283
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:437
Foam::VectorSpace< Vector< scalar >, scalar, 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:411