1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011-2011 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 "cylinderAnnulusToCell.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 defineTypeNameAndDebug(cylinderAnnulusToCell, 0);
35 addToRunTimeSelectionTable(topoSetSource, cylinderAnnulusToCell, word);
36 addToRunTimeSelectionTable(topoSetSource, cylinderAnnulusToCell, istream);
40 Foam::topoSetSource::addToUsageTable Foam::cylinderAnnulusToCell::usage_
42 cylinderAnnulusToCell::typeName,
43 "\n Usage: cylinderAnnulusToCell (p1X p1Y p1Z) (p2X p2Y p2Z)"
44 " outerRadius innerRadius\n\n"
45 " Select all cells with cell centre within bounding cylinder annulus\n\n"
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 void Foam::cylinderAnnulusToCell::combine(topoSet& set, const bool add) const
53 const vector axis = p2_ - p1_;
54 const scalar orad2 = sqr(outerRadius_);
55 const scalar irad2 = sqr(innerRadius_);
56 const scalar magAxis2 = magSqr(axis);
58 const pointField& ctrs = mesh_.cellCentres();
62 vector d = ctrs[cellI] - p1_;
63 scalar magD = d & axis;
65 if ((magD > 0) && (magD < magAxis2))
67 scalar d2 = (d & d) - sqr(magD)/magAxis2;
68 if ((d2 < orad2) && (d2 > irad2))
70 addOrDelete(set, cellI, add);
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 Foam::cylinderAnnulusToCell::cylinderAnnulusToCell
84 const scalar outerRadius,
85 const scalar innerRadius
91 outerRadius_(outerRadius),
92 innerRadius_(innerRadius)
96 Foam::cylinderAnnulusToCell::cylinderAnnulusToCell
99 const dictionary& dict
103 p1_(dict.lookup("p1")),
104 p2_(dict.lookup("p2")),
105 outerRadius_(readScalar(dict.lookup("outerRadius"))),
106 innerRadius_(readScalar(dict.lookup("innerRadius")))
110 // Construct from Istream
111 Foam::cylinderAnnulusToCell::cylinderAnnulusToCell
113 const polyMesh& mesh,
120 outerRadius_(readScalar(checkIs(is))),
121 innerRadius_(readScalar(checkIs(is)))
125 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
127 Foam::cylinderAnnulusToCell::~cylinderAnnulusToCell()
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 void Foam::cylinderAnnulusToCell::applyToSet
135 const topoSetSource::setAction action,
139 if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
141 Info<< " Adding cells with centre within cylinder annulus,"
143 << p1_ << ", p2 = " << p2_ << " and outer radius = " << outerRadius_
144 << " and inner radius = " << innerRadius_
149 else if (action == topoSetSource::DELETE)
151 Info<< " Removing cells with centre within cylinder, with p1 = "
152 << p1_ << ", p2 = " << p2_ << " and outer radius = " << outerRadius_
153 << " and inner radius " << innerRadius_
161 // ************************************************************************* //