sampledCuttingSurface.H
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) 2018-2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::sampledCuttingSurface
28
29Description
30 A surface define by using an input surface to cut the mesh cells
31
32Usage
33 Dictionary controls:
34 \table
35 Property | Description | Required | Default
36 type | surfaceCut | yes |
37 surfaceType | type of surface | yes |
38 surfaceName | name of surface in \c triSurface/ | no | dict name
39 triangulate | triangulate faces | no | true
40 bounds | limit with bounding box | no |
41 zone | limit to cell zone (name or regex) | no |
42 zones | limit to cell zones (names, regexs) | no |
43 \endtable
44
45Note
46 No attempt at resolving degenerate cases.
47 Since the cut faces can be quite ugly, they will often be triangulated.
48
49SourceFiles
50 sampledCuttingSurface.C
51 sampledCuttingSurfaceTemplates.C
52
53\*---------------------------------------------------------------------------*/
54
55#ifndef sampledCuttingSurface_H
56#define sampledCuttingSurface_H
57
58#include "sampledSurface.H"
59#include "cuttingSurface.H"
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63namespace Foam
64{
65
66/*---------------------------------------------------------------------------*\
67 Class sampledCuttingSurface Declaration
68\*---------------------------------------------------------------------------*/
69
70class sampledCuttingSurface
71:
72 public sampledSurface,
73 public cuttingSurface
74{
75 // Private Data
76
77 //- The zone or zones in which cutting is to occur
78 wordRes zoneNames_;
79
80 //- Optional bounding box to trim against
81 const boundBox bounds_;
82
83 //- Triangulate faces or not
84 const bool triangulate_;
85
86 //- Track if the surface needs an update
87 mutable bool needsUpdate_;
88
89
90 // Private Member Functions
91
92 //- Define cell selection from zones and bounding box.
93 // Optionally check and warn if the bounding box
94 // does not overlap with the mesh (or submesh)
95 bitSet cellSelection(const bool warn=false) const;
96
97
98 //- Sample volume field onto surface faces
99 template<class Type>
100 tmp<Field<Type>> sampleOnFaces
101 (
102 const interpolation<Type>& sampler
103 ) const;
104
105 //- Interpolate volume field onto surface points
106 template<class Type>
107 tmp<Field<Type>> sampleOnPoints
108 (
109 const interpolation<Type>& interpolator
110 ) const;
111
112
113public:
114
115 //- Runtime type information
116 TypeName("surfaceCut");
117
118
119 // Constructors
120
121 //- Construct from components
123 (
124 const polyMesh& mesh,
125 const word& surfaceType,
126 const word& surfaceName,
127 const bool triangulate = true,
128 const boundBox& bounds = boundBox::invertedBox
129 );
130
131 //- Construct from dictionary
133 (
134 const word& defaultSurfaceName,
135 const polyMesh& mesh,
136 const dictionary& dict
137 );
138
139
140 //- Destructor
141 virtual ~sampledCuttingSurface() = default;
142
143
144 // Member Functions
145
146 //- Does the surface need an update?
147 virtual bool needsUpdate() const;
148
149 //- Mark the surface as needing an update.
150 // May also free up unneeded data.
151 // Return false if surface was already marked as expired.
152 virtual bool expire();
153
154 //- Update the surface as required.
155 // Do nothing (and return false) if no update was needed
156 virtual bool update();
157
158
159 //- Points of surface
160 virtual const pointField& points() const
161 {
162 return cuttingSurface::points();
163 }
164
165 //- Faces of surface
166 virtual const faceList& faces() const
167 {
169 }
170
171 //- Per-face zone/region information
172 // Could instead return meshCells or cellZoneId of the meshCells.
173 virtual const labelList& zoneIds() const
174 {
175 return labelList::null();
176 }
177
178 //- Face area magnitudes
179 virtual const vectorField& Sf() const
181 return cuttingSurface::Sf();
182 }
183
184 //- Face area magnitudes
185 virtual const scalarField& magSf() const
186 {
187 return cuttingSurface::magSf();
188 }
189
190 //- Face centres
191 virtual const vectorField& Cf() const
192 {
193 return cuttingSurface::Cf();
194 }
195
196 //- For each face, the original cell in mesh
198
200 // Sample
201
202 //- Sample volume field onto surface faces
204 (
206 ) const;
207
208 //- Sample volume field onto surface faces
210 (
211 const interpolation<vector>& sampler
212 ) const;
213
214 //- Sample volume field onto surface faces
216 (
217 const interpolation<sphericalTensor>& sampler
218 ) const;
219
220 //- Sample volume field onto surface faces
222 (
223 const interpolation<symmTensor>& sampler
224 ) const;
225
226 //- Sample volume field onto surface faces
228 (
229 const interpolation<tensor>& sampler
230 ) const;
231
232
233 // Interpolate
234
235 //- Interpolate volume field onto surface points
237 (
238 const interpolation<scalar>& interpolator
239 ) const;
240
241 //- Interpolate volume field onto surface points
243 (
244 const interpolation<vector>& interpolator
245 ) const;
246
247 //- Interpolate volume field onto surface points
249 (
250 const interpolation<sphericalTensor>& interpolator
251 ) const;
252
253 //- Interpolate volume field onto surface points
255 (
256 const interpolation<symmTensor>& interpolator
257 ) const;
258
259 //- Interpolate volume field onto surface points
261 (
262 const interpolation<tensor>& interpolator
263 ) const;
264};
265
266
267// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268
269} // End namespace Foam
270
271// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272
273#ifdef NoRepository
275#endif
276
277// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278
279#endif
280
281// ************************************************************************* //
Minimal example by using system/controlDict.functions:
static const List< label > & null()
Return a null List.
Definition: ListI.H:109
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
const List< Face > & surfFaces() const
Return const access to the faces.
surfaceTopo surfaceType() const
Calculate surface type formed by patch.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Constructs a cutting surface through a mesh.
const word & surfaceName() const
The name of the underlying searchableSurface.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A surface define by using an input surface to cut the mesh cells.
virtual const pointField & points() const
Points of surface.
TypeName("surfaceCut")
Runtime type information.
virtual const faceList & faces() const
Faces of surface.
virtual const vectorField & Cf() const
Face centres.
virtual const labelList & zoneIds() const
Per-face zone/region information.
virtual const scalarField & magSf() const
Face area magnitudes.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool update()
Update the surface as required.
virtual const vectorField & Sf() const
Face area magnitudes.
virtual ~sampledCuttingSurface()=default
Destructor.
const labelList & meshCells() const
For each face, the original cell in mesh.
An abstract class for surfaces with sampling.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
bool interpolate() const noexcept
Same as isPointData()
A class for managing temporary objects.
Definition: tmp.H:65
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
const pointField & points
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
Field< vector > vectorField
Specialisation of Field<T> for vector.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73