sampledPlane.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2022 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
27Class
28 Foam::sampledPlane
29
30Description
31 A sampledSurface defined by a plane which \a cuts the mesh using the
32 cuttingPlane alorithm. The surface is triangulated by default.
33
34 This is often embedded as part of a sampled surfaces function object.
35
36Usage
37 Example of function object partial specification:
38 \verbatim
39 surfaces
40 {
41 surface1
42 {
43 type plane;
44 point ...;
45 normal ...;
46 }
47 }
48 \endverbatim
49
50 Where the sub-entries comprise:
51 \table
52 Property | Description | Reqd | Default
53 type | plane | yes |
54 planeType | Plane description (pointAndNormal etc) | no |
55 triangulate | triangulate faces | no | true
56 bounds | limit with bounding box | no |
57 zone | limit to cell zone (name or regex) | no |
58 zones | limit to cell zones (names, regexs) | no |
59 coordinateSystem | define plane within given cartesian system | no |
60 transform | define plane within given cartesian system | no |
61 \endtable
62
63SeeAlso
64 Foam::coordinateSystem
65 Foam::cuttingPlane
66 Foam::plane
67
68Note
69 Does not actually cut until update() called.
70 The keyword \c zones has priority over \c zone.
71
72SourceFiles
73 sampledPlane.C
74
75\*---------------------------------------------------------------------------*/
76
77#ifndef Foam_sampledPlane_H
78#define Foam_sampledPlane_H
79
80#include "sampledSurface.H"
81#include "cuttingPlane.H"
82
83// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84
85namespace Foam
86{
87
88/*---------------------------------------------------------------------------*\
89 Class sampledPlane Declaration
90\*---------------------------------------------------------------------------*/
91
92class sampledPlane
93:
94 public sampledSurface,
95 public cuttingPlane
96{
97 // Private data
98
99 //- The zone or zones in which cutting is to occur
100 wordRes zoneNames_;
101
102 //- Optional bounding box to trim against
103 const boundBox bounds_;
104
105 //- Triangulate faces or not
106 const bool triangulate_;
107
108 //- Track if the surface needs an update
109 mutable bool needsUpdate_;
110
111
112 // Private Member Functions
113
114 //- Define plane from dictionary entry and apply any coordinate
115 //- transformations
116 static plane definePlane(const polyMesh& mesh, const dictionary& dict);
117
118 //- Define cell selection from zones and bounding box.
119 // Optionally check and warn if the plane does not intersect
120 // with the bounds of the mesh (or submesh) or if the bounding box
121 // does not overlap with the mesh (or submesh)
122 bitSet cellSelection(const bool warn=false) const;
123
124
125 //- Sample volume field onto surface faces
126 template<class Type>
127 tmp<Field<Type>> sampleOnFaces
128 (
129 const interpolation<Type>& sampler
130 ) const;
131
132 //- Interpolate volume field onto surface points
133 template<class Type>
134 tmp<Field<Type>> sampleOnPoints
135 (
136 const interpolation<Type>& interpolator
137 ) const;
138
139
140public:
141
142 //- Runtime type information
143 TypeName("sampledPlane");
144
145
146 // Constructors
147
148 //- Construct from components
150 (
151 const word& name,
152 const polyMesh& mesh,
153 const plane& planeDesc,
154 const wordRes& zones = wordRes(),
155 const bool triangulate = true
156 );
157
158 //- Construct from dictionary
160 (
161 const word& name,
162 const polyMesh& mesh,
163 const dictionary& dict
164 );
165
166
167 //- Destructor
168 virtual ~sampledPlane() = default;
169
170
171 // Member Functions
172
173 //- Does the surface need an update?
174 virtual bool needsUpdate() const;
175
176 //- Mark the surface as needing an update.
177 // May also free up unneeded data.
178 // Return false if surface was already marked as expired.
179 virtual bool expire();
180
181 //- Update the surface as required.
182 // Do nothing (and return false) if no update was needed
183 virtual bool update();
184
185
186 //- Points of surface
187 virtual const pointField& points() const
188 {
189 return cuttingPlane::points();
190 }
191
192 //- Faces of surface
193 virtual const faceList& faces() const
194 {
196 }
197
198 //- Per-face zone/region information
199 // Could instead return meshCells or cellZoneId of the meshCells.
200 virtual const labelList& zoneIds() const
201 {
202 return labelList::null();
203 }
204
205 //- Face area magnitudes
206 virtual const vectorField& Sf() const
207 {
208 return cuttingPlane::Sf();
209 }
210
211 //- Face area magnitudes
212 virtual const scalarField& magSf() const
213 {
214 return cuttingPlane::magSf();
215 }
216
217 //- Face centres
218 virtual const vectorField& Cf() const
219 {
220 return cuttingPlane::Cf();
221 }
222
223 //- For each face, the original cell in mesh
225
226
227 // Sample
228
229 //- Sample volume field onto surface faces
230 virtual tmp<scalarField> sample
232 const interpolation<scalar>& sampler
233 ) const;
234
235 //- Sample volume field onto surface faces
238 const interpolation<vector>& sampler
239 ) const;
240
241 //- Sample volume field onto surface faces
243 (
245 ) const;
246
247 //- Sample volume field onto surface faces
249 (
251 ) const;
252
253 //- Sample volume field onto surface faces
255 (
257 ) const;
258
259
260 // Interpolate
261
262 //- Interpolate volume field onto surface points
264 (
265 const interpolation<scalar>& interpolator
266 ) const;
267
268 //- Interpolate volume field onto surface points
270 (
271 const interpolation<vector>& interpolator
272 ) const;
273
274 //- Interpolate volume field onto surface points
276 (
277 const interpolation<sphericalTensor>& interpolator
278 ) const;
279
280 //- Interpolate volume field onto surface points
282 (
283 const interpolation<symmTensor>& interpolator
284 ) const;
285
286 //- Interpolate volume field onto surface points
288 (
289 const interpolation<tensor>& interpolator
290 ) const;
291
292
293 // Output
294
295 //- Print information
296 virtual void print(Ostream& os, int level=0) const;
297};
298
299
300// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301
302} // End namespace Foam
303
304// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305
306#ifdef NoRepository
307 #include "sampledPlaneTemplates.C"
308#endif
309
310// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311
312#endif
313
314// ************************************************************************* //
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
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
Constructs cutting plane through a mesh.
Definition: cuttingPlane.H:62
const plane & planeDesc() const
The plane used.
Definition: cuttingPlane.H:172
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
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:95
plane()
Construct zero-initialised.
Definition: planeI.H:30
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const labelList & meshCells() const
For each face, the original cell in mesh.
A sampledSurface defined by a plane which cuts the mesh using the cuttingPlane alorithm....
Definition: sampledPlane.H:140
virtual const pointField & points() const
Points of surface.
Definition: sampledPlane.H:231
TypeName("sampledPlane")
Runtime type information.
virtual void print(Ostream &os, int level=0) const
Print information.
Definition: sampledPlane.C:338
virtual const faceList & faces() const
Faces of surface.
Definition: sampledPlane.H:237
virtual const vectorField & Cf() const
Face centres.
Definition: sampledPlane.H:262
virtual const labelList & zoneIds() const
Per-face zone/region information.
Definition: sampledPlane.H:244
virtual const scalarField & magSf() const
Face area magnitudes.
Definition: sampledPlane.H:256
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPlane.C:212
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPlane.C:206
virtual bool update()
Update the surface as required.
Definition: sampledPlane.C:227
virtual const vectorField & Sf() const
Face area magnitudes.
Definition: sampledPlane.H:250
virtual ~sampledPlane()=default
Destructor.
An abstract class for surfaces with sampling.
const word & name() const noexcept
Name of surface.
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
OBJstream os(runTime.globalPath()/outputName)
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