1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "searchablePlate.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "SortableList.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 defineTypeNameAndDebug(searchablePlate, 0);
37 addToRunTimeSelectionTable(searchableSurface, searchablePlate, dict);
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 Foam::direction Foam::searchablePlate::calcNormal(const point& span)
46 direction normalDir = 3;
48 for (direction dir = 0; dir < vector::nComponents; dir++)
52 FatalErrorIn("searchablePlate::calcNormal()")
53 << "Span should have two positive and one zero entry. Now:"
54 << span << exit(FatalError);
56 else if (span[dir] < VSMALL)
64 // Multiple zero entries. Flag and exit.
73 FatalErrorIn("searchablePlate::calcNormal()")
74 << "Span should have one and only zero entry. Now:" << span
82 // Returns miss or hit with face (always 0)
83 Foam::pointIndexHit Foam::searchablePlate::findNearest
86 const scalar nearestDistSqr
89 // For every component direction can be
90 // left of min, right of max or inbetween.
91 // - outside points: project first one x plane (either min().x()
92 // or max().x()), then onto y plane and finally z. You should be left
93 // with intersection point
94 // - inside point: find nearest side (compare to mid point). Project onto
97 // Project point on plane.
98 pointIndexHit info(true, sample, 0);
99 info.rawPoint()[normalDir_] = origin_[normalDir_];
101 // Clip to edges if outside
102 for (direction dir = 0; dir < vector::nComponents; dir++)
104 if (dir != normalDir_)
106 if (info.rawPoint()[dir] < origin_[dir])
108 info.rawPoint()[dir] = origin_[dir];
110 else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
112 info.rawPoint()[dir] = origin_[dir]+span_[dir];
117 // Check if outside. Optimisation: could do some checks on distance already
118 // on components above
119 if (magSqr(info.rawPoint() - sample) > nearestDistSqr)
129 Foam::pointIndexHit Foam::searchablePlate::findLine
142 const vector dir(end-start);
144 if (mag(dir[normalDir_]) < VSMALL)
151 scalar t = (origin_[normalDir_]-start[normalDir_]) / dir[normalDir_];
160 info.rawPoint() = start+t*dir;
161 info.rawPoint()[normalDir_] = origin_[normalDir_];
164 for (direction dir = 0; dir < vector::nComponents; dir++)
166 if (dir != normalDir_)
168 if (info.rawPoint()[dir] < origin_[dir])
174 else if (info.rawPoint()[dir] > origin_[dir]+span_[dir])
188 treeBoundBox bb(origin_, origin_+span_);
189 bb.min()[normalDir_] -= 1E-6;
190 bb.max()[normalDir_] += 1E-6;
192 if (!bb.contains(info.hitPoint()))
194 FatalErrorIn("searchablePlate::findLine(..)")
195 << "bb:" << bb << endl
196 << "origin_:" << origin_ << endl
197 << "span_:" << span_ << endl
198 << "normalDir_:" << normalDir_ << endl
199 << "hitPoint:" << info.hitPoint()
200 << abort(FatalError);
208 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
210 Foam::searchablePlate::searchablePlate
217 searchableSurface(io),
220 normalDir_(calcNormal(span_))
224 Info<< "searchablePlate::searchablePlate :"
225 << " origin:" << origin_
226 << " origin+span:" << origin_+span_
227 << " normal:" << vector::componentNames[normalDir_]
233 Foam::searchablePlate::searchablePlate
236 const dictionary& dict
239 searchableSurface(io),
240 origin_(dict.lookup("origin")),
241 span_(dict.lookup("span")),
242 normalDir_(calcNormal(span_))
246 Info<< "searchablePlate::searchablePlate :"
247 << " origin:" << origin_
248 << " origin+span:" << origin_+span_
249 << " normal:" << vector::componentNames[normalDir_]
255 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
257 Foam::searchablePlate::~searchablePlate()
261 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
263 const Foam::wordList& Foam::searchablePlate::regions() const
265 if (regions_.empty())
268 regions_[0] = "region0";
274 void Foam::searchablePlate::findNearest
276 const pointField& samples,
277 const scalarField& nearestDistSqr,
278 List<pointIndexHit>& info
281 info.setSize(samples.size());
285 info[i] = findNearest(samples[i], nearestDistSqr[i]);
290 void Foam::searchablePlate::findLine
292 const pointField& start,
293 const pointField& end,
294 List<pointIndexHit>& info
297 info.setSize(start.size());
301 info[i] = findLine(start[i], end[i]);
306 void Foam::searchablePlate::findLineAny
308 const pointField& start,
309 const pointField& end,
310 List<pointIndexHit>& info
313 findLine(start, end, info);
317 void Foam::searchablePlate::findLineAll
319 const pointField& start,
320 const pointField& end,
321 List<List<pointIndexHit> >& info
324 List<pointIndexHit> nearestInfo;
325 findLine(start, end, nearestInfo);
327 info.setSize(start.size());
330 if (nearestInfo[pointI].hit())
332 info[pointI].setSize(1);
333 info[pointI][0] = nearestInfo[pointI];
337 info[pointI].clear();
343 void Foam::searchablePlate::getRegion
345 const List<pointIndexHit>& info,
349 region.setSize(info.size());
354 void Foam::searchablePlate::getNormal
356 const List<pointIndexHit>& info,
360 normal.setSize(info.size());
361 normal = vector::zero;
364 normal[i][normalDir_] = 1.0;
369 void Foam::searchablePlate::getVolumeType
371 const pointField& points,
372 List<volumeType>& volType
377 "searchableCollection::getVolumeType(const pointField&"
378 ", List<volumeType>&) const"
379 ) << "Volume type not supported for plate."
384 // ************************************************************************* //