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-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 "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 {
39 }
40 
41 
42 // Angle for polys to be considered splitHexes.
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::shapeToCell::combine(topoSet& set, const bool add) const
49 {
50  if (type_ == "splitHex")
51  {
52  for (label celli = 0; celli < mesh_.nCells(); celli++)
53  {
54  cellFeatures superCell(mesh_, featureCos, celli);
55 
56  if (hexMatcher().isA(superCell.faces()))
57  {
58  addOrDelete(set, celli, add);
59  }
60  }
61  }
62  else
63  {
64  const cellModel& wantedModel = *(cellModeller::lookup(type_));
65 
67 
68  forAll(cellShapes, celli)
69  {
70  if (cellShapes[celli].model() == wantedModel)
71  {
72  addOrDelete(set, celli, add);
73  }
74  }
75  }
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
82 (
83  const polyMesh& mesh,
84  const word& type
85 )
86 :
87  topoSetSource(mesh),
88  type_(type)
89 {
90  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
91  {
93  << "Illegal cell type " << type_ << exit(FatalError);
94  }
95 }
96 
97 
99 (
100  const polyMesh& mesh,
101  const dictionary& dict
102 )
103 :
104  topoSetSource(mesh),
105  type_(dict.lookup("type"))
106 {
107  if (!cellModeller::lookup(type_) && (type_ != "splitHex"))
108  {
110  << "Illegal cell type " << type_ << exit(FatalError);
111  }
112 }
113 
114 
115 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
116 
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 (
125  const topoSetSource::setAction action,
126  topoSet& set
127 ) const
128 {
129  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
130  {
131  Info<< " Adding all cells of type " << type_ << " ..." << endl;
132 
133  combine(set, true);
134  }
135  else if (action == topoSetSource::DELETE)
136  {
137  Info<< " Removing all cells of type " << type_ << " ..." << endl;
138 
139  combine(set, false);
140  }
141 }
142 
143 
144 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
label nCells() const
const cellShapeList & cellShapes() const
Return cell shapes.
A topoSetSource to select cells based on cell shape.
Definition: shapeToCell.H:54
static scalar featureCos
Cos of feature angle for polyHedral to be splitHex.
Definition: shapeToCell.H:77
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
Definition: shapeToCell.C:124
virtual ~shapeToCell()
Destructor.
Definition: shapeToCell.C:117
shapeToCell(const polyMesh &mesh, const word &type)
Construct from components.
Definition: shapeToCell.C:82
Base class of a source for a topoSet.
Definition: topoSetSource.H:64
void addOrDelete(topoSet &set, const label celli, const bool) const
Add (if bool) celli to set or delete celli from set.
Definition: topoSetSource.C:96
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:83
const polyMesh & mesh_
Definition: topoSetSource.H:99
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:65
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const cellShapeList & cellShapes
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:43
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:139
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
treeBoundBox combine(const treeBoundBox &a, const treeBoundBox &b)
Definition: patchToPatch.C:78
dimensionedScalar cos(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dictionary dict
Unit conversion functions.