shapeToCell.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-2018 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 "shapeToCell.H"
27 #include "polyMesh.H"
28 #include "unitConversion.H"
29 #include "hexMatcher.H"
30 #include "cellFeatures.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(shapeToCell, 0);
38  addToRunTimeSelectionTable(topoSetSource, shapeToCell, word);
39  addToRunTimeSelectionTable(topoSetSource, shapeToCell, istream);
40 }
41 
42 
43 Foam::topoSetSource::addToUsageTable Foam::shapeToCell::usage_
44 (
45  shapeToCell::typeName,
46  "\n Usage: shapeToCell tet|pyr|prism|hex|tetWedge|wedge|splitHex\n\n"
47  " Select all cells of given cellShape.\n"
48  " (splitHex hardcoded with internal angle < 10 degrees)\n"
49 );
50 
51 
52 // Angle for polys to be considered splitHexes.
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
58 void Foam::shapeToCell::combine(topoSet& set, const bool add) const
59 {
60  if (type_ == "splitHex")
61  {
62  for (label celli = 0; celli < mesh_.nCells(); celli++)
63  {
64  cellFeatures superCell(mesh_, featureCos, celli);
65 
66  if (hexMatcher().isA(superCell.faces()))
67  {
68  addOrDelete(set, celli, add);
69  }
70  }
71  }
72  else
73  {
74  const cellModel& wantedModel = *(cellModeller::lookup(type_));
75 
76  const cellShapeList& cellShapes = mesh_.cellShapes();
77 
78  forAll(cellShapes, celli)
79  {
80  if (cellShapes[celli].model() == wantedModel)
81  {
82  addOrDelete(set, celli, add);
83  }
84  }
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const polyMesh& mesh,
94  const word& type
95 )
96 :
97  topoSetSource(mesh),
98  type_(type)
99 {
100  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
101  {
103  << "Illegal cell type " << type_ << exit(FatalError);
104  }
105 }
106 
107 
109 (
110  const polyMesh& mesh,
111  const dictionary& dict
112 )
113 :
114  topoSetSource(mesh),
115  type_(dict.lookup("type"))
116 {
117  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
118  {
120  << "Illegal cell type " << type_ << exit(FatalError);
121  }
122 }
123 
124 
126 (
127  const polyMesh& mesh,
128  Istream& is
129 )
130 :
131  topoSetSource(mesh),
132  type_(checkIs(is))
133 {
134  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
135  {
137  << "Illegal cell type " << type_ << exit(FatalError);
138  }
139 }
140 
141 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
142 
144 {}
145 
146 
147 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148 
150 (
151  const topoSetSource::setAction action,
152  topoSet& set
153 ) const
154 {
155  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
156  {
157  Info<< " Adding all cells of type " << type_ << " ..." << endl;
158 
159  combine(set, true);
160  }
161  else if (action == topoSetSource::DELETE)
162  {
163  Info<< " Removing all cells of type " << type_ << " ..." << endl;
164 
165  combine(set, false);
166  }
167 }
168 
169 
170 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or nullptr.
Definition: cellModeller.C:100
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Unit conversion functions.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Definition: ListListOps.C:34
Macros for easy insertion into run-time selection tables.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Base class of a source for a topoSet.
Definition: topoSetSource.H:63
virtual ~shapeToCell()
Destructor.
Definition: shapeToCell.C:143
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:43
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:82
static scalar featureCos
Cos of feature angle for polyHedral to be splitHex.
Definition: shapeToCell.H:79
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
Definition: shapeToCell.C:150
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
shapeToCell(const polyMesh &mesh, const word &type)
Construct from components.
Definition: shapeToCell.C:92
Class with constructor to add usage string to table.
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576