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 | Copyright (C) 2011-2016 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"
31 
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 defineTypeNameAndDebug(shapeToCell, 0);
40 
41 addToRunTimeSelectionTable(topoSetSource, shapeToCell, word);
42 
43 addToRunTimeSelectionTable(topoSetSource, shapeToCell, istream);
44 
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 (type_ == "splitHex")
66  {
67  for (label celli = 0; celli < mesh_.nCells(); celli++)
68  {
69  cellFeatures superCell(mesh_, featureCos, celli);
70 
71  if (hexMatcher().isA(superCell.faces()))
72  {
73  addOrDelete(set, celli, add);
74  }
75  }
76  }
77  else
78  {
79  const cellModel& wantedModel = *(cellModeller::lookup(type_));
80 
81  const cellShapeList& cellShapes = mesh_.cellShapes();
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 
96 // Construct from components
98 (
99  const polyMesh& mesh,
100  const word& type
101 )
102 :
103  topoSetSource(mesh),
104  type_(type)
105 {
106  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
107  {
109  << "Illegal cell type " << type_ << exit(FatalError);
110  }
111 }
112 
113 
114 // Construct from dictionary
116 (
117  const polyMesh& mesh,
118  const dictionary& dict
119 )
120 :
121  topoSetSource(mesh),
122  type_(dict.lookup("type"))
123 {
124  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
125  {
127  << "Illegal cell type " << type_ << exit(FatalError);
128  }
129 }
130 
131 
132 // Construct from Istream
134 (
135  const polyMesh& mesh,
136  Istream& is
137 )
138 :
139  topoSetSource(mesh),
140  type_(checkIs(is))
141 {
142  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
143  {
145  << "Illegal cell type " << type_ << exit(FatalError);
146  }
147 }
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
152 {}
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
158 (
159  const topoSetSource::setAction action,
160  topoSet& set
161 ) const
162 {
163  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
164  {
165  Info<< " Adding all cells of type " << type_ << " ..." << endl;
166 
167  combine(set, true);
168  }
169  else if (action == topoSetSource::DELETE)
170  {
171  Info<< " Removing all cells of type " << type_ << " ..." << endl;
172 
173  combine(set, false);
174  }
175 }
176 
177 
178 // ************************************************************************* //
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:88
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:253
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:151
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:158
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:98
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