1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2010 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
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 "cylinderToCell.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 defineTypeNameAndDebug(cylinderToCell, 0);
35 addToRunTimeSelectionTable(topoSetSource, cylinderToCell, word);
36 addToRunTimeSelectionTable(topoSetSource, cylinderToCell, istream);
40 Foam::topoSetSource::addToUsageTable Foam::cylinderToCell::usage_
42 cylinderToCell::typeName,
43 "\n Usage: cylinderToCell (p1X p1Y p1Z) (p2X p2Y p2Z) radius\n\n"
44 " Select all cells with cell centre within bounding cylinder\n\n"
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 void Foam::cylinderToCell::combine(topoSet& set, const bool add) const
52 const vector axis = p2_ - p1_;
53 const scalar rad2 = sqr(radius_);
54 const scalar magAxis2 = magSqr(axis);
56 const pointField& ctrs = mesh_.cellCentres();
60 vector d = ctrs[cellI] - p1_;
61 scalar magD = d & axis;
63 if ((magD > 0) && (magD < magAxis2))
65 scalar d2 = (d & d) - sqr(magD)/magAxis2;
68 addOrDelete(set, cellI, add);
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 Foam::cylinderToCell::cylinderToCell
92 Foam::cylinderToCell::cylinderToCell
95 const dictionary& dict
99 p1_(dict.lookup("p1")),
100 p2_(dict.lookup("p2")),
101 radius_(readScalar(dict.lookup("radius")))
105 // Construct from Istream
106 Foam::cylinderToCell::cylinderToCell
108 const polyMesh& mesh,
115 radius_(readScalar(checkIs(is)))
119 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
121 Foam::cylinderToCell::~cylinderToCell()
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 void Foam::cylinderToCell::applyToSet
129 const topoSetSource::setAction action,
133 if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
135 Info<< " Adding cells with centre within cylinder, with p1 = "
136 << p1_ << ", p2 = " << p2_ << " and radius = " << radius_ << endl;
140 else if (action == topoSetSource::DELETE)
142 Info<< " Removing cells with centre within cylinder, with p1 = "
143 << p1_ << ", p2 = " << p2_ << " and radius = " << radius_ << endl;
150 // ************************************************************************* //