cylinderAnnulusToCell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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
19  for more details.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "cylinderAnnulusToCell.H"
27 #include "polyMesh.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(cylinderAnnulusToCell, 0);
35  addToRunTimeSelectionTable(topoSetSource, cylinderAnnulusToCell, word);
36 }
37 
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 void Foam::cylinderAnnulusToCell::combine(topoSet& set, const bool add) const
42 {
43  const vector axis = point2_ - point1_;
44  const scalar orad2 = sqr(outerRadius_);
45  const scalar irad2 = sqr(innerRadius_);
46  const scalar magAxis2 = magSqr(axis);
47 
48  const pointField& ctrs = mesh_.cellCentres();
49 
50  forAll(ctrs, celli)
51  {
52  vector d = ctrs[celli] - point1_;
53  scalar magD = d & axis;
54 
55  if ((magD > 0) && (magD < magAxis2))
56  {
57  scalar d2 = (d & d) - sqr(magD)/magAxis2;
58  if ((d2 < orad2) && (d2 > irad2))
59  {
60  addOrDelete(set, celli, add);
61  }
62  }
63  }
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
70 (
71  const polyMesh& mesh,
72  const vector& point1,
73  const vector& point2,
74  const scalar outerRadius,
75  const scalar innerRadius
76 )
77 :
78  topoSetSource(mesh),
79  point1_(point1),
80  point2_(point2),
81  outerRadius_(outerRadius),
82  innerRadius_(innerRadius)
83 {}
84 
85 
87 (
88  const polyMesh& mesh,
89  const dictionary& dict
90 )
91 :
92  topoSetSource(mesh),
93  point1_(dict.lookupBackwardsCompatible<point>({"point1", "p1"})),
94  point2_(dict.lookupBackwardsCompatible<point>({"point2", "p2"})),
95  outerRadius_(dict.lookup<scalar>("outerRadius")),
96  innerRadius_(dict.lookup<scalar>("innerRadius"))
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
101 
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 (
110  const topoSetSource::setAction action,
111  topoSet& set
112 ) const
113 {
114  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
115  {
116  Info<< " Adding cells with centre within cylinder annulus,"
117  << " with point1 = "
118  << point1_ << ", point2 = " << point2_ << " and outer radius = "
119  << outerRadius_ << " and inner radius = " << innerRadius_ << endl;
120 
121  combine(set, true);
122  }
123  else if (action == topoSetSource::DELETE)
124  {
125  Info<< " Removing cells with centre within cylinder, with point1 = "
126  << point1_ << ", point2 = " << point2_ << " and outer radius = "
127  << outerRadius_ << " and inner radius " << innerRadius_ << endl;
128 
129  combine(set, false);
130  }
131 }
132 
133 
134 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Macros for easy insertion into run-time selection tables.
Base class of a source for a topoSet.
Definition: topoSetSource.H:63
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
cylinderAnnulusToCell(const polyMesh &mesh, const vector &point1, const vector &point2, const scalar outerRadius, const scalar innerRadius)
Construct from components.
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:82
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
dimensioned< scalar > magSqr(const dimensioned< Type > &)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:875
Namespace for OpenFOAM.
treeBoundBox combine(const treeBoundBox &a, const treeBoundBox &b)
Definition: patchToPatch.C:77
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864