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-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Class
28  Foam::sampledPlane
29 
30 Description
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 
36 Usage
37  Example of function object partial specification:
38  \verbatim
39  surfaces
40  {
41  surface1
42  {
43  type plane;
44  planeType pointAndNormal;
45  pointAndNormalDict
46  {
47  ...
48  }
49  }
50  }
51  \endverbatim
52 
53  Where the sub-entries comprise:
54  \table
55  Property | Description | Required | Default
56  type | plane | yes |
57  planeType | plane description (pointAndNormal etc) | yes |
58  triangulate | triangulate faces | no | true
59  bounds | limit with bounding box | no |
60  zone | limit to cell zone (name or regex) | no |
61  zones | limit to cell zones (names, regexs) | no |
62  coordinateSystem | define plane within given coordinate system | no |
63  \endtable
64 
65 SeeAlso
66  Foam::coordinateSystem
67  Foam::cuttingPlane
68  Foam::plane
69 
70 Note
71  Does not actually cut until update() called.
72  The keyword \c zones has priority over \c zone.
73 
74 SourceFiles
75  sampledPlane.C
76 
77 \*---------------------------------------------------------------------------*/
78 
79 #ifndef sampledPlane_H
80 #define sampledPlane_H
81 
82 #include "sampledSurface.H"
83 #include "cuttingPlane.H"
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 namespace Foam
88 {
89 
90 /*---------------------------------------------------------------------------*\
91  Class sampledPlane Declaration
92 \*---------------------------------------------------------------------------*/
93 
94 class sampledPlane
95 :
96  public sampledSurface,
97  public cuttingPlane
98 {
99  // Private data
100 
101  //- The zone or zones in which cutting is to occur
102  wordRes zoneNames_;
103 
104  //- Optional bounding box to trim against
105  const boundBox bounds_;
106 
107  //- Triangulate faces or not
108  const bool triangulate_;
109 
110  //- Track if the surface needs an update
111  mutable bool needsUpdate_;
112 
113 
114  // Private Member Functions
115 
116  //- Define cell selection from zones and bounding box.
117  // Optionally check and warn if the plane does not intersect
118  // with the bounds of the mesh (or submesh) or if the bounding box
119  // does not overlap with the mesh (or submesh)
120  bitSet cellSelection(const bool warn=false) const;
121 
122 
123  //- Sample volume field onto surface faces
124  template<class Type>
125  tmp<Field<Type>> sampleOnFaces
126  (
127  const interpolation<Type>& sampler
128  ) const;
129 
130  //- Interpolate volume field onto surface points
131  template<class Type>
132  tmp<Field<Type>> sampleOnPoints
133  (
134  const interpolation<Type>& interpolator
135  ) const;
136 
137 
138 public:
139 
140  //- Runtime type information
141  TypeName("sampledPlane");
142 
143 
144  // Constructors
145 
146  //- Construct from components
148  (
149  const word& name,
150  const polyMesh& mesh,
151  const plane& planeDesc,
152  const wordRes& zones = wordRes(),
153  const bool triangulate = true
154  );
155 
156  //- Construct from dictionary
158  (
159  const word& name,
160  const polyMesh& mesh,
161  const dictionary& dict
162  );
163 
164 
165  //- Destructor
166  virtual ~sampledPlane() = default;
167 
168 
169  // Member Functions
170 
171  //- Does the surface need an update?
172  virtual bool needsUpdate() const;
173 
174  //- Mark the surface as needing an update.
175  // May also free up unneeded data.
176  // Return false if surface was already marked as expired.
177  virtual bool expire();
178 
179  //- Update the surface as required.
180  // Do nothing (and return false) if no update was needed
181  virtual bool update();
182 
183 
184  //- Points of surface
185  virtual const pointField& points() const
186  {
187  return cuttingPlane::points();
188  }
189 
190  //- Faces of surface
191  virtual const faceList& faces() const
192  {
193  return cuttingPlane::surfFaces();
194  }
195 
196  //- Per-face zone/region information
197  // Could instead return meshCells or cellZoneId of the meshCells.
198  virtual const labelList& zoneIds() const
199  {
200  return labelList::null();
201  }
202 
203  //- Face area magnitudes
204  virtual const vectorField& Sf() const
205  {
206  return cuttingPlane::Sf();
207  }
208 
209  //- Face area magnitudes
210  virtual const scalarField& magSf() const
211  {
212  return cuttingPlane::magSf();
213  }
214 
215  //- Face centres
216  virtual const vectorField& Cf() const
217  {
218  return cuttingPlane::Cf();
219  }
220 
221  //- For each face, the original cell in mesh
223 
224 
225  // Sample
226 
227  //- Sample volume field onto surface faces
228  virtual tmp<scalarField> sample
229  (
230  const interpolation<scalar>& sampler
231  ) const;
232 
233  //- Sample volume field onto surface faces
234  virtual tmp<vectorField> sample
235  (
236  const interpolation<vector>& sampler
237  ) const;
238 
239  //- Sample volume field onto surface faces
241  (
242  const interpolation<sphericalTensor>& sampler
243  ) const;
244 
245  //- Sample volume field onto surface faces
247  (
248  const interpolation<symmTensor>& sampler
249  ) const;
250 
251  //- Sample volume field onto surface faces
252  virtual tmp<tensorField> sample
253  (
254  const interpolation<tensor>& sampler
255  ) const;
256 
257 
258  // Interpolate
259 
260  //- Interpolate volume field onto surface points
262  (
263  const interpolation<scalar>& interpolator
264  ) const;
265 
266  //- Interpolate volume field onto surface points
268  (
269  const interpolation<vector>& interpolator
270  ) const;
271 
272  //- Interpolate volume field onto surface points
274  (
275  const interpolation<sphericalTensor>& interpolator
276  ) const;
277 
278  //- Interpolate volume field onto surface points
280  (
281  const interpolation<symmTensor>& interpolator
282  ) const;
283 
284  //- Interpolate volume field onto surface points
286  (
287  const interpolation<tensor>& interpolator
288  ) const;
289 
290 
291  // Output
292 
293  //- Print information
294  virtual void print(Ostream& os, int level=0) const;
295 };
296 
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 } // End namespace Foam
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 #ifdef NoRepository
305  #include "sampledPlaneTemplates.C"
306 #endif
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 #endif
311 
312 // ************************************************************************* //
Foam::sampledPlane::sampledPlane
sampledPlane(const word &name, const polyMesh &mesh, const plane &planeDesc, const wordRes &zones=wordRes(), const bool triangulate=true)
Construct from components.
Definition: sampledPlane.C:69
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::sampledPlane::update
virtual bool update()
Update the surface as required.
Definition: sampledPlane.C:193
Foam::cuttingSurfaceBase::meshCells
const labelList & meshCells() const
The mesh cells cut.
Definition: cuttingSurfaceBase.H:197
Foam::List::null
static const List< T > & null()
Return a null List.
Definition: ListI.H:109
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::sampledPlane::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPlane.C:178
Foam::MeshedSurface::triangulate
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
Definition: MeshedSurface.C:1015
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::cuttingPlane
Constructs cutting plane through a mesh.
Definition: cuttingPlane.H:58
Foam::MeshedSurface::Sf
const vectorField & Sf() const
Face area vectors (normals)
Definition: MeshedSurface.H:435
Foam::sampledPlane::magSf
virtual const scalarField & magSf() const
Face area magnitudes.
Definition: sampledPlane.H:249
Foam::MeshedSurface::surfFaces
const List< Face > & surfFaces() const
Return const access to the faces.
Definition: MeshedSurface.H:413
Foam::sampledPlane::sample
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
Definition: sampledPlane.C:216
Foam::sampledPlane::TypeName
TypeName("sampledPlane")
Runtime type information.
cuttingPlane.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::sampledPlane::faces
virtual const faceList & faces() const
Faces of surface.
Definition: sampledPlane.H:230
Foam::sampledPlane::points
virtual const pointField & points() const
Points of surface.
Definition: sampledPlane.H:224
Foam::sampledPlane::print
virtual void print(Ostream &os, int level=0) const
Print information.
Definition: sampledPlane.C:304
Foam::cuttingPlane::planeDesc
const plane & planeDesc() const
The plane used.
Definition: cuttingPlane.H:172
Foam::MeshedSurface::magSf
const scalarField & magSf() const
Face area magnitudes.
Definition: MeshedSurface.H:441
Foam::sampledPlane::~sampledPlane
virtual ~sampledPlane()=default
Destructor.
Foam::sampledSurface::interpolate
bool interpolate() const noexcept
Same as isPointData()
Definition: sampledSurface.H:598
Foam::sampledPlane::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPlane.C:172
Foam::Field< vector >
Foam::sampledPlane::Cf
virtual const vectorField & Cf() const
Face centres.
Definition: sampledPlane.H:255
sampledSurface.H
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:121
Foam::sampledPlane::zoneIds
virtual const labelList & zoneIds() const
Per-face zone/region information.
Definition: sampledPlane.H:237
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::sampledSurface::name
const word & name() const noexcept
Name of surface.
Definition: sampledSurface.H:322
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sampledPlane::Sf
virtual const vectorField & Sf() const
Face area magnitudes.
Definition: sampledPlane.H:243
Foam::List< face >
Foam::sampledSurface::mesh
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
Definition: sampledSurface.H:316
sampledPlaneTemplates.C
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::MeshedSurface::Cf
const vectorField & Cf() const
Face centres.
Definition: MeshedSurface.H:447
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::sampledPlane
A sampledSurface defined by a plane which cuts the mesh using the cuttingPlane alorithm....
Definition: sampledPlane.H:133