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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::sampledPlane
26 
27 Description
28  A sampledSurface defined by a plane which 'cuts' the mesh using the
29  cuttingPlane alorithm. The plane is triangulated by default.
30 
31 Note
32  Does not actually cut until update() called.
33 
34 SourceFiles
35  sampledPlane.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef sampledPlane_H
40 #define sampledPlane_H
41 
42 #include "sampledSurface.H"
43 #include "cuttingPlane.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class sampledPlane Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 class sampledPlane
55 :
56  public sampledSurface,
57  public cuttingPlane
58 {
59  // Private data
60 
61  //- If restricted to zones, name of this zone or a regular expression
62  keyType zoneKey_;
63 
64  //- Triangulated faces or keep faces as is
65  const bool triangulate_;
66 
67  //- Track if the surface needs an update
68  mutable bool needsUpdate_;
69 
70  // Private Member Functions
71 
72  //- Sample field on faces
73  template<class Type>
74  tmp<Field<Type>> sampleField
75  (
77  ) const;
78 
79 
80  template<class Type>
82  interpolateField(const interpolation<Type>&) const;
83 
84 
85 public:
86 
87  //- Runtime type information
88  TypeName("sampledPlane");
89 
90 
91  // Constructors
92 
93  //- Construct from components
95  (
96  const word& name,
97  const polyMesh& mesh,
98  const plane& planeDesc,
99  const keyType& zoneKey = word::null,
100  const bool triangulate = true
101  );
102 
103  //- Construct from dictionary
105  (
106  const word& name,
107  const polyMesh& mesh,
108  const dictionary& dict
109  );
110 
111 
112  //- Destructor
113  virtual ~sampledPlane();
114 
115 
116  // Member Functions
117 
118  //- Does the surface need an update?
119  virtual bool needsUpdate() const;
120 
121  //- Mark the surface as needing an update.
122  // May also free up unneeded data.
123  // Return false if surface was already marked as expired.
124  virtual bool expire();
125 
126  //- Update the surface as required.
127  // Do nothing (and return false) if no update was needed
128  virtual bool update();
129 
130 
131  //- Points of surface
132  virtual const pointField& points() const
133  {
134  return cuttingPlane::points();
135  }
136 
137  //- Faces of surface
138  virtual const faceList& faces() const
139  {
140  return cuttingPlane::faces();
141  }
142 
143  //- For every face original cell in mesh
144  const labelList& meshCells() const
145  {
146  return cuttingPlane::cutCells();
147  }
148 
149  //- Sample field on surface
150  virtual tmp<scalarField> sample
151  (
152  const volScalarField&
153  ) const;
154 
155 
156  //- Sample field on surface
157  virtual tmp<vectorField> sample
158  (
159  const volVectorField&
160  ) const;
161 
162  //- Sample field on surface
164  (
166  ) const;
167 
168  //- Sample field on surface
170  (
171  const volSymmTensorField&
172  ) const;
173 
174  //- Sample field on surface
175  virtual tmp<tensorField> sample
176  (
177  const volTensorField&
178  ) const;
179 
180 
181  //- Interpolate field on surface
183  (
184  const interpolation<scalar>&
185  ) const;
186 
187 
188  //- Interpolate field on surface
190  (
191  const interpolation<vector>&
192  ) const;
193 
194  //- Interpolate field on surface
196  (
198  ) const;
199 
200  //- Interpolate field on surface
202  (
204  ) const;
205 
206  //- Interpolate field on surface
208  (
209  const interpolation<tensor>&
210  ) const;
211 
212  //- Write
213  virtual void print(Ostream&) const;
214 };
215 
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 } // End namespace Foam
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #ifdef NoRepository
224  #include "sampledPlaneTemplates.C"
225 #endif
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 #endif
230 
231 // ************************************************************************* //
A class for handling keywords in dictionaries.
Definition: keyType.H:64
virtual void print(Ostream &) const
Write.
Definition: sampledPlane.C:251
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An abstract class for surfaces with sampling.
const word & name() const
Name of surface.
bool interpolate() const
Interpolation requested for surface.
Generic GeometricField class.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPlane.C:116
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledPlane.C:163
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPlane.C:110
Constructs plane through mesh.
Definition: cuttingPlane.H:60
A class for handling words, derived from string.
Definition: word.H:59
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
const Field< PointType > & points() const
Return reference to global points.
static const word null
An empty word.
Definition: word.H:77
virtual bool update()
Update the surface as required.
Definition: sampledPlane.C:131
const polyMesh & mesh() const
Access to the underlying mesh.
virtual const pointField & points() const
Points of surface.
Definition: sampledPlane.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const List< face > & faces() const
Return const access to the faces.
A sampledSurface defined by a plane which &#39;cuts&#39; the mesh using the cuttingPlane alorithm. The plane is triangulated by default.
Definition: sampledPlane.H:53
virtual ~sampledPlane()
Destructor.
Definition: sampledPlane.C:104
Abstract base class for interpolation.
const labelList & cutCells() const
Return List of cells cut by the plane.
Definition: cuttingPlane.H:159
const plane & planeDesc() const
Return plane used.
Definition: cuttingPlane.H:153
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
TypeName("sampledPlane")
Runtime type information.
const labelList & meshCells() const
For every face original cell in mesh.
Definition: sampledPlane.H:143
virtual const faceList & faces() const
Faces of surface.
Definition: sampledPlane.H:137
Namespace for OpenFOAM.
sampledPlane(const word &name, const polyMesh &mesh, const plane &planeDesc, const keyType &zoneKey=word::null, const bool triangulate=true)
Construct from components.
Definition: sampledPlane.C:44