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"
33#include "cellFeatures.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
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
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")
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:
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// ************************************************************************* //
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:64
reference ref() const
A reference to the entry (Error if not found)
Definition: dictionary.H:225
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
static bool test(const UList< face > &faces)
Definition: hexMatcher.C:87
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const cellShapeList & cellShapes() const
Return cell shapes.
label nCells() const noexcept
Number of mesh cells.
A topoSetCellSource to select cells based on the type of their cell shapes.
Definition: shapeToCell.H:159
static scalar featureCos
Cos of feature angle for polyHedral to be splitHex.
Definition: shapeToCell.H:183
virtual void applyToSet(const topoSetSource::setAction action, topoSet &set) const
Apply specified action to the topoSet.
Definition: shapeToCell.C:143
The topoSetCellSource is a intermediate class for handling topoSet sources for selecting cells.
Class with constructor to add usage string to table.
Base class of a source for a topoSet.
Definition: topoSetSource.H:68
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.
const polyMesh & mesh_
Reference to the mesh.
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:67
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
cellShapeList cellShapes
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Namespace for OpenFOAM.
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:45
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionedScalar cos(const dimensionedScalar &ds)
dict add("bounds", meshBb)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Unit conversion functions.