Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38namespace Foam
39{
45}
46
47
48Foam::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.
58Foam::scalar Foam::shapeToCell::featureCos = Foam::cos(degToRad(10.0));
59
60
61// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62
63void 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
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
97(
98 const polyMesh& mesh,
99 const word& shapeName
100)
101:
103 shape_(shapeName)
104{
105 if (!cellModel::ptr(shape_) && shape_ != "splitHex")
108 << "Illegal cell shape " << shape_ << exit(FatalError);
109 }
110}
111
112
114(
115 const polyMesh& mesh,
116 const dictionary& dict
118:
120 shape_(dict.getCompat<word>("shape", {{"type", 1806}}))
121{}
122
123
125(
126 const polyMesh& mesh,
127 Istream& is
128)
129:
130 topoSetCellSource(mesh),
131 shape_(checkIs(is))
132{
133 if (!cellModel::ptr(shape_) && shape_ != "splitHex")
134 {
136 << "Illegal cell shape " << shape_ << exit(FatalError);
137 }
138}
139
140
141// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142
144(
145 const topoSetSource::setAction action,
146 topoSet& set
147) const
148{
149 if (action == topoSetSource::ADD || action == topoSetSource::NEW)
150 {
151 if (verbose_)
152 {
153 Info<< " Adding all " << shape_ << " cells ..." << endl;
154 }
155
156 combine(set, true);
157 }
158 else if (action == topoSetSource::SUBTRACT)
159 {
160 if (verbose_)
161 {
162 Info<< " Removing all " << shape_ << " cells ..." << endl;
163 }
164
165 combine(set, false);
166 }
167}
168
169
170// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
static const cellModel * ptr(const modelType model)
Look up pointer to cellModel by enumeration, or nullptr on failure.
Definition cellModels.C:113
static const cellModel & ref(const modelType model)
Look up reference to cellModel by enumeration. Fatal on failure.
Definition cellModels.C:150
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for HEX. (6 quad).
Definition hexMatcher.C:80
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
label nCells() const noexcept
Number of mesh cells.
A topoSetCellSource to select cells based on the type of their cell shapes.
static scalar featureCos
Cos of feature angle for polyHedral to be splitHex.
shapeToCell(const polyMesh &mesh, const word &shapeName)
Construct from components.
Definition shapeToCell.C:90
virtual void applyToSet(const topoSetSource::setAction action, topoSet &set) const
Apply specified action to the topoSet.
The topoSetCellSource is a intermediate class for handling topoSet sources for selecting cells.
topoSetCellSource(const polyMesh &mesh)
Construct from mesh.
Class with constructor to add usage string to table.
Base class of a source for a topoSet.
void addOrDelete(topoSet &set, const label id, const bool add) const
Add or delete id from set. Add when 'add' is true.
setAction
Enumeration defining various actions.
@ SUBTRACT
Subtract elements from current set.
@ ADD
Add elements to current set.
@ NEW
Create a new set and ADD elements to it.
bool verbose_
Output verbosity (default: true).
const polyMesh & mesh() const noexcept
Reference to the mesh.
const polyMesh & mesh_
Reference to the mesh.
static Istream & checkIs(Istream &is)
Check state of stream.
General set of labels of mesh quantity (points, cells, faces).
Definition topoSet.H:63
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
cellShapeList cellShapes
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
List< cellShape > cellShapeList
List of cellShape.
dict add("bounds", meshBb)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.