shapeToCell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "shapeToCell.H"
30 #include "polyMesh.H"
31 #include "unitConversion.H"
32 #include "hexMatcher.H"
33 #include "cellFeatures.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(shapeToCell, 0);
41  addToRunTimeSelectionTable(topoSetSource, shapeToCell, word);
42  addToRunTimeSelectionTable(topoSetSource, shapeToCell, istream);
43  addToRunTimeSelectionTable(topoSetCellSource, shapeToCell, word);
44  addToRunTimeSelectionTable(topoSetCellSource, shapeToCell, istream);
45 }
46 
47 
48 Foam::topoSetSource::addToUsageTable Foam::shapeToCell::usage_
49 (
50  shapeToCell::typeName,
51  "\n Usage: shapeToCell tet|pyr|prism|hex|tetWedge|wedge|splitHex\n\n"
52  " Select all cells of given cellShape.\n"
53  " (splitHex hardcoded with internal angle < 10 degrees)\n"
54 );
55 
56 
57 // Angle for polys to be considered splitHexes.
59 
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62 
63 void Foam::shapeToCell::combine(topoSet& set, const bool add) const
64 {
65  if (shape_ == "splitHex")
66  {
67  for (label celli = 0; celli < mesh_.nCells(); ++celli)
68  {
69  cellFeatures superCell(mesh_, featureCos, celli);
70 
71  if (hexMatcher::test(superCell.faces()))
72  {
73  addOrDelete(set, celli, add);
74  }
75  }
76  }
77  else
78  {
79  const cellModel& wantedModel = cellModel::ref(shape_);
80 
82 
83  forAll(cellShapes, celli)
84  {
85  if (cellShapes[celli].model() == wantedModel)
86  {
87  addOrDelete(set, celli, add);
88  }
89  }
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const polyMesh& mesh,
99  const word& shapeName
100 )
101 :
102  topoSetCellSource(mesh),
103  shape_(shapeName)
104 {
105  if (!cellModel::ptr(shape_) && shape_ != "splitHex")
106  {
108  << "Illegal cell shape " << shape_ << exit(FatalError);
109  }
110 }
111 
112 
114 (
115  const polyMesh& mesh,
116  const dictionary& dict
117 )
118 :
119  shapeToCell(mesh, dict.getCompat<word>("shape", {{"type", 1806}}))
120 {}
121 
122 
124 (
125  const polyMesh& mesh,
126  Istream& is
127 )
128 :
129  topoSetCellSource(mesh),
130  shape_(checkIs(is))
131 {
132  if (!cellModel::ptr(shape_) && shape_ != "splitHex")
133  {
135  << "Illegal cell shape " << shape_ << exit(FatalError);
136  }
137 }
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 (
144  const topoSetSource::setAction action,
145  topoSet& set
146 ) const
147 {
148  if (action == topoSetSource::ADD || action == topoSetSource::NEW)
149  {
150  if (verbose_)
151  {
152  Info<< " Adding all " << shape_ << " cells ..." << endl;
153  }
154 
155  combine(set, true);
156  }
157  else if (action == topoSetSource::SUBTRACT)
158  {
159  if (verbose_)
160  {
161  Info<< " Removing all " << shape_ << " cells ..." << endl;
162  }
163 
164  combine(set, false);
165  }
166 }
167 
168 
169 // ************************************************************************* //
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Create a new set and ADD elements to it.
Add elements to current set.
Unit conversion functions.
const cellShapeList & cellShapes() const
Return cell shapes.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
List< cellShape > cellShapeList
List of cellShape.
Definition: cellShapeList.H:32
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for HEX. (6 quad)
Definition: hexMatcher.C:80
void addOrDelete(topoSet &set, const label id, const bool add) const
Add or delete id from set. Add when &#39;add&#39; is true.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static const cellModel & ref(const modelType model)
Look up reference to cellModel by enumeration. Fatal on failure.
Definition: cellModels.C:150
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Definition: ListListOps.C:62
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
setAction
Enumeration defining various actions.
static const cellModel * ptr(const modelType model)
Look up pointer to cellModel by enumeration, or nullptr on failure.
Definition: cellModels.C:113
cellShapeList cellShapes
shapeToCell(const polyMesh &mesh, const word &shapeName)
Construct from components.
Definition: shapeToCell.C:90
const polyMesh & mesh_
Reference to the mesh.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
static scalar featureCos
Cos of feature angle for polyHedral to be splitHex.
Definition: shapeToCell.H:187
defineTypeNameAndDebug(combustionModel, 0)
A topoSetCellSource to select cells based on the type of their cell shapes.
Definition: shapeToCell.H:152
Subtract elements from current set.
Class with constructor to add usage string to table.
label nCells() const noexcept
Number of mesh cells.
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
virtual void applyToSet(const topoSetSource::setAction action, topoSet &set) const
Apply specified action to the topoSet.
Definition: shapeToCell.C:136
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)