1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "searchablePlate.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "SortableList.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 defineTypeNameAndDebug(searchablePlate, 0);
36 addToRunTimeSelectionTable(searchableSurface, searchablePlate, dict);
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 Foam::direction Foam::searchablePlate::calcNormal(const point& span)
45 direction normalDir = 3;
47 for (direction dir = 0; dir < vector::nComponents; dir++)
51 FatalErrorIn("searchablePlate::calcNormal()")
52 << "Span should have two positive and one zero entry. Now:"
53 << span << exit(FatalError);
55 else if (span[dir] < VSMALL)
63 // Multiple zero entries. Flag and exit.
72 FatalErrorIn("searchablePlate::calcNormal()")
73 << "Span should have one and only zero entry. Now:" << span
81 // Returns miss or hit with face (always 0)
82 Foam::pointIndexHit Foam::searchablePlate::findNearest
85 const scalar nearestDistSqr
88 // For every component direction can be
89 // left of min, right of max or inbetween.
90 // - outside points: project first one x plane (either min().x()
91 // or max().x()), then onto y plane and finally z. You should be left
92 // with intersection point
93 // - inside point: find nearest side (compare to mid point). Project onto
96 // Project point on plane.
97 pointIndexHit info(true, sample, 0);
98 info.rawPoint()[normalDir_] = origin_[normalDir_];
100 // Clip to edges if outside
101 for (direction dir = 0; dir < vector::nComponents; dir++)
103 if (dir != normalDir_)
105 if (info.rawPoint()[dir] < origin_[dir])
107 info.rawPoint()[dir] = origin_[dir];
109 else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
111 info.rawPoint()[dir] = origin_[dir]+span_[dir];
116 // Check if outside. Optimisation: could do some checks on distance already
117 // on components above
118 if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
128 Foam::pointIndexHit Foam::searchablePlate::findLine
141 const vector dir(end-start);
143 if (mag(dir[normalDir_]) < VSMALL)
150 scalar t = (origin_[normalDir_]-start[normalDir_]) / dir[normalDir_];
159 info.rawPoint() = start+t*dir;
160 info.rawPoint()[normalDir_] = origin_[normalDir_];
163 for (direction dir = 0; dir < vector::nComponents; dir++)
165 if (dir != normalDir_)
167 if (info.rawPoint()[dir] < origin_[dir])
173 else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
187 treeBoundBox bb(origin_, origin_+span_);
188 bb.min()[normalDir_] -= 1E-6;
189 bb.max()[normalDir_] += 1E-6;
191 if (!bb.contains(info.hitPoint()))
193 FatalErrorIn("searchablePlate::findLine(..)")
194 << "bb:" << bb << endl
195 << "origin_:" << origin_ << endl
196 << "span_:" << span_ << endl
197 << "normalDir_:" << normalDir_ << endl
198 << "hitPoint:" << info.hitPoint()
199 << abort(FatalError);
207 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
209 Foam::searchablePlate::searchablePlate
216 searchableSurface(io),
219 normalDir_(calcNormal(span_))
223 Info<< "searchablePlate::searchablePlate :"
224 << " origin:" << origin_
225 << " origin+span:" << origin_+span_
226 << " normal:" << vector::componentNames[normalDir_]
232 Foam::searchablePlate::searchablePlate
235 const dictionary& dict
238 searchableSurface(io),
239 origin_(dict.lookup("origin")),
240 span_(dict.lookup("span")),
241 normalDir_(calcNormal(span_))
245 Info<< "searchablePlate::searchablePlate :"
246 << " origin:" << origin_
247 << " origin+span:" << origin_+span_
248 << " normal:" << vector::componentNames[normalDir_]
254 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
256 Foam::searchablePlate::~searchablePlate()
260 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
262 const Foam::wordList& Foam::searchablePlate::regions() const
264 if (regions_.empty())
267 regions_[0] = "region0";
273 void Foam::searchablePlate::findNearest
275 const pointField& samples,
276 const scalarField& nearestDistSqr,
277 List<pointIndexHit>& info
280 info.setSize(samples.size());
284 info[i] = findNearest(samples[i], nearestDistSqr[i]);
289 void Foam::searchablePlate::findLine
291 const pointField& start,
292 const pointField& end,
293 List<pointIndexHit>& info
296 info.setSize(start.size());
300 info[i] = findLine(start[i], end[i]);
305 void Foam::searchablePlate::findLineAny
307 const pointField& start,
308 const pointField& end,
309 List<pointIndexHit>& info
312 findLine(start, end, info);
316 void Foam::searchablePlate::findLineAll
318 const pointField& start,
319 const pointField& end,
320 List<List<pointIndexHit> >& info
323 List<pointIndexHit> nearestInfo;
324 findLine(start, end, nearestInfo);
326 info.setSize(start.size());
329 if (nearestInfo[pointI].hit())
331 info[pointI].setSize(1);
332 info[pointI][0] = nearestInfo[pointI];
336 info[pointI].clear();
342 void Foam::searchablePlate::getRegion
344 const List<pointIndexHit>& info,
348 region.setSize(info.size());
353 void Foam::searchablePlate::getNormal
355 const List<pointIndexHit>& info,
359 normal.setSize(info.size());
360 normal = vector::zero;
363 normal[i][normalDir_] = 1.0;
368 void Foam::searchablePlate::getVolumeType
370 const pointField& points,
371 List<volumeType>& volType
376 "searchableCollection::getVolumeType(const pointField&"
377 ", List<volumeType>&) const"
378 ) << "Volume type not supported for plate."
383 // ************************************************************************* //