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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "searchablePlate.H"
30 #include "SortableList.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(searchablePlate, 0);
38  (
39  searchableSurface,
40  searchablePlate,
41  dict
42  );
44  (
45  searchableSurface,
46  searchablePlate,
47  dict,
48  plate
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 Foam::direction Foam::searchablePlate::calcNormal(const point& span)
56 {
57  direction normalDir = 3;
58 
59  for (direction dir = 0; dir < vector::nComponents; ++dir)
60  {
61  if (span[dir] < 0)
62  {
64  << "Span should have two positive and one zero entry. Now:"
65  << span << exit(FatalError);
66  }
67  else if (span[dir] < VSMALL)
68  {
69  if (normalDir == 3)
70  {
71  normalDir = dir;
72  }
73  else
74  {
75  // Multiple zero entries. Flag and exit.
76  normalDir = 3;
77  break;
78  }
79  }
80  }
81 
82  if (normalDir == 3)
83  {
85  << "Span should have two positive and one zero entry. Now:"
86  << span << 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 {
233  if (debug)
234  {
236  << " origin:" << origin_
237  << " origin+span:" << origin_+span_
238  << " normal:" << vector::componentNames[normalDir_]
239  << endl;
240  }
241 
242  bounds() = boundBox(origin_, origin_ + span_);
243 }
244 
245 
246 Foam::searchablePlate::searchablePlate
247 (
248  const IOobject& io,
249  const dictionary& dict
250 )
251 :
252  searchableSurface(io),
253  origin_(dict.get<point>("origin")),
254  span_(dict.get<vector>("span")),
255  normalDir_(calcNormal(span_))
256 {
257  if (debug)
258  {
260  << " origin:" << origin_
261  << " origin+span:" << origin_+span_
262  << " normal:" << vector::componentNames[normalDir_]
263  << endl;
264  }
265 
266  bounds() = boundBox(origin_, origin_ + span_);
267 }
268 
269 
270 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
271 
273 {
274  if (regions_.empty())
275  {
276  regions_.setSize(1);
277  regions_[0] = "region0";
278  }
279  return regions_;
280 }
281 
282 
284 {
285  return tmp<pointField>::New(1, origin_ + 0.5*span_);
286 }
287 
288 
290 (
291  pointField& centres,
292  scalarField& radiusSqr
293 ) const
294 {
295  centres.setSize(1);
296  centres[0] = origin_ + 0.5*span_;
297 
298  radiusSqr.setSize(1);
299  radiusSqr[0] = Foam::magSqr(0.5*span_);
300 
301  // Add a bit to make sure all points are tested inside
302  radiusSqr += Foam::sqr(SMALL);
303 }
304 
305 
307 {
308  auto tpts = tmp<pointField>::New(4);
309  auto& pts = tpts.ref();
310 
311  pts[0] = origin_;
312  pts[2] = origin_ + span_;
313 
314  if (span_.x() < span_.y() && span_.x() < span_.z())
315  {
316  pts[1] = origin_ + point(0, span_.y(), 0);
317  pts[3] = origin_ + point(0, 0, span_.z());
318  }
319  else if (span_.y() < span_.z())
320  {
321  pts[1] = origin_ + point(span_.x(), 0, 0);
322  pts[3] = origin_ + point(0, 0, span_.z());
323  }
324  else
325  {
326  pts[1] = origin_ + point(span_.x(), 0, 0);
327  pts[3] = origin_ + point(0, span_.y(), 0);
328  }
329 
330  return tpts;
331 }
332 
333 
335 {
336  return
337  (
338  (origin_.x() + span_.x()) >= bb.min().x()
339  && origin_.x() <= bb.max().x()
340  && (origin_.y() + span_.y()) >= bb.min().y()
341  && origin_.y() <= bb.max().y()
342  && (origin_.z() + span_.z()) >= bb.min().z()
343  && origin_.z() <= bb.max().z()
344  );
345 }
346 
347 
348 void Foam::searchablePlate::findNearest
349 (
350  const pointField& samples,
351  const scalarField& nearestDistSqr,
352  List<pointIndexHit>& info
353 ) const
354 {
355  info.setSize(samples.size());
356 
357  forAll(samples, i)
358  {
359  info[i] = findNearest(samples[i], nearestDistSqr[i]);
360  }
361 }
362 
363 
364 void Foam::searchablePlate::findLine
365 (
366  const pointField& start,
367  const pointField& end,
368  List<pointIndexHit>& info
369 ) const
370 {
371  info.setSize(start.size());
372 
373  forAll(start, i)
374  {
375  info[i] = findLine(start[i], end[i]);
376  }
377 }
378 
379 
381 (
382  const pointField& start,
383  const pointField& end,
384  List<pointIndexHit>& info
385 ) const
386 {
387  findLine(start, end, info);
388 }
389 
390 
392 (
393  const pointField& start,
394  const pointField& end,
396 ) const
397 {
398  List<pointIndexHit> nearestInfo;
399  findLine(start, end, nearestInfo);
400 
401  info.setSize(start.size());
402  forAll(info, pointi)
403  {
404  if (nearestInfo[pointi].hit())
405  {
406  info[pointi].setSize(1);
407  info[pointi][0] = nearestInfo[pointi];
408  }
409  else
410  {
411  info[pointi].clear();
412  }
413  }
414 }
415 
416 
418 (
419  const List<pointIndexHit>& info,
420  labelList& region
421 ) const
422 {
423  region.setSize(info.size());
424  region = 0;
425 }
426 
427 
429 (
430  const List<pointIndexHit>& info,
431  vectorField& normal
432 ) const
433 {
434  normal.setSize(info.size());
435  normal = Zero;
436  forAll(normal, i)
437  {
438  normal[i][normalDir_] = 1.0;
439  }
440 }
441 
442 
444 (
445  const pointField& points,
446  List<volumeType>& volType
447 ) const
448 {
450  << "Volume type not supported for plate."
451  << exit(FatalError);
452 }
453 
454 
455 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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:78
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
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:392
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
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:334
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:90
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:290
Foam::searchablePlate::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: searchablePlate.C:306
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 >
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:115
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:84
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
Foam::Vector< scalar >
Foam::List< word >
Foam::searchablePlate::coordinates
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
Definition: searchablePlate.C:283
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:381
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
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:272
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:429
Foam::searchablePlate::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
Definition: searchablePlate.C:290
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:444
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:418