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 cuttingPlane. Triangulated by default.
29 
30 Note
31  Does not actually cut until update() called.
32 
33 SourceFiles
34  sampledPlane.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef sampledPlane_H
39 #define sampledPlane_H
40 
41 #include "sampledSurface.H"
42 #include "cuttingPlane.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class sampledPlane Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 class sampledPlane
54 :
55  public sampledSurface,
56  public cuttingPlane
57 {
58  // Private data
59 
60  //- If restricted to zones, name of this zone or a regular expression
61  keyType zoneKey_;
62 
63  //- Triangulated faces or keep faces as is
64  const bool triangulate_;
65 
66  //- Track if the surface needs an update
67  mutable bool needsUpdate_;
68 
69  // Private Member Functions
70 
71  //- Sample field on faces
72  template<class Type>
73  tmp<Field<Type>> sampleField
74  (
76  ) const;
77 
78 
79  template<class Type>
81  interpolateField(const interpolation<Type>&) const;
82 
83 
84 public:
85 
86  //- Runtime type information
87  TypeName("sampledPlane");
88 
89 
90  // Constructors
91 
92  //- Construct from components
94  (
95  const word& name,
96  const polyMesh& mesh,
97  const plane& planeDesc,
98  const keyType& zoneKey = word::null,
99  const bool triangulate = true
100  );
101 
102  //- Construct from dictionary
104  (
105  const word& name,
106  const polyMesh& mesh,
107  const dictionary& dict
108  );
109 
110 
111  //- Destructor
112  virtual ~sampledPlane();
113 
114 
115  // Member Functions
116 
117  //- Does the surface need an update?
118  virtual bool needsUpdate() const;
119 
120  //- Mark the surface as needing an update.
121  // May also free up unneeded data.
122  // Return false if surface was already marked as expired.
123  virtual bool expire();
124 
125  //- Update the surface as required.
126  // Do nothing (and return false) if no update was needed
127  virtual bool update();
128 
129 
130  //- Points of surface
131  virtual const pointField& points() const
132  {
133  return cuttingPlane::points();
134  }
135 
136  //- Faces of surface
137  virtual const faceList& faces() const
138  {
139  return cuttingPlane::faces();
140  }
141 
142  //- For every face original cell in mesh
143  const labelList& meshCells() const
144  {
145  return cuttingPlane::cutCells();
146  }
147 
148  //- Sample field on surface
149  virtual tmp<scalarField> sample
150  (
151  const volScalarField&
152  ) const;
153 
154 
155  //- Sample field on surface
156  virtual tmp<vectorField> sample
157  (
158  const volVectorField&
159  ) const;
160 
161  //- Sample field on surface
163  (
165  ) const;
166 
167  //- Sample field on surface
169  (
170  const volSymmTensorField&
171  ) const;
172 
173  //- Sample field on surface
174  virtual tmp<tensorField> sample
175  (
176  const volTensorField&
177  ) const;
178 
179 
180  //- Interpolate field on surface
182  (
183  const interpolation<scalar>&
184  ) const;
185 
186 
187  //- Interpolate field on surface
189  (
190  const interpolation<vector>&
191  ) const;
192 
193  //- Interpolate field on surface
195  (
197  ) const;
198 
199  //- Interpolate field on surface
201  (
203  ) const;
204 
205  //- Interpolate field on surface
207  (
208  const interpolation<tensor>&
209  ) const;
210 
211  //- Write
212  virtual void print(Ostream&) const;
213 };
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 #ifdef NoRepository
223  #include "sampledPlaneTemplates.C"
224 #endif
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
A class for handling keywords in dictionaries.
Definition: keyType.H:64
virtual const faceList & faces() const
Faces of surface.
Definition: sampledPlane.H:136
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const labelList & cutCells() const
Return List of cells cut by the plane.
Definition: cuttingPlane.H:158
An abstract class for surfaces with sampling.
virtual const pointField & points() const
Points of surface.
Definition: sampledPlane.H:130
const List< face > & faces() const
Return const access to the faces.
virtual void print(Ostream &) const
Write.
Definition: sampledPlane.C:251
const plane & planeDesc() const
Return plane used.
Definition: cuttingPlane.H:152
const Field< PointType > & points() const
Return reference to global points.
const labelList & meshCells() const
For every face original cell in mesh.
Definition: sampledPlane.H:142
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
const polyMesh & mesh() const
Access to the underlying mesh.
Constructs plane through mesh.
Definition: cuttingPlane.H:60
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPlane.C:110
A class for handling words, derived from string.
Definition: word.H:59
const word & name() const
Name of surface.
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
static const word null
An empty word.
Definition: word.H:77
virtual bool update()
Update the surface as required.
Definition: sampledPlane.C:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A sampledSurface defined by a cuttingPlane. Triangulated by default.
Definition: sampledPlane.H:52
bool interpolate() const
Interpolation requested for surface.
virtual ~sampledPlane()
Destructor.
Definition: sampledPlane.C:104
Abstract base class for interpolation.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:54
TypeName("sampledPlane")
Runtime type information.
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledPlane.C:163
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