sampledPlane.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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::sampledSurfaces::plane
26 
27 Description
28  A sampledSurface defined by a plane which 'cuts' the mesh using the
29  cuttingPlane algorithm. The plane is triangulated by default.
30 
31  Example:
32  \verbatim
33  {
34  type plane;
35  planeType pointAndNormal;
36  point (0 0 0);
37  normal (0 0 1);
38  interpolate yes;
39  }
40  \endverbatim
41 
42 Usage
43  \table
44  Property | Description | Required | Default value
45  planeType | the method of specification of the plane | yes |
46  triangulate | triangulate the output | no | yes
47  interpolate | interpolate values to the surface points | no | no
48  \endtable
49 
50 See also
51  Foam::plane
52 
53 SourceFiles
54  sampledPlane.C
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #ifndef sampledPlane_H
59 #define sampledPlane_H
60 
61 #include "sampledSurface.H"
62 #include "cuttingPlane.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 namespace sampledSurfaces
69 {
70 
71 /*---------------------------------------------------------------------------*\
72  Class plane Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class plane
76 :
77  public sampledSurface,
78  public cuttingPlane
79 {
80  // Private Data
81 
82  //- If restricted to zones, name of this zone or a regular expression
83  wordRe zoneKey_;
84 
85  //- Triangulated faces or keep faces as is
86  const bool triangulate_;
87 
88  //- Track if the surface needs an update
89  mutable bool needsUpdate_;
90 
91  // Private Member Functions
92 
93  //- Sample field on faces
94  template<class Type>
95  tmp<Field<Type>> sampleField
96  (
98  ) const;
99 
100 
101  template<class Type>
103  interpolateField(const interpolation<Type>&) const;
104 
105 
106 public:
107 
108  //- Runtime type information
109  TypeName("plane");
110 
111 
112  // Constructors
113 
114  //- Construct from dictionary
115  plane
116  (
117  const word& name,
118  const polyMesh& mesh,
119  const dictionary& dict
120  );
121 
122 
123  //- Destructor
124  virtual ~plane();
125 
126 
127  // Member Functions
128 
129  //- Does the surface need an update?
130  virtual bool needsUpdate() const;
131 
132  //- Mark the surface as needing an update.
133  // May also free up unneeded data.
134  // Return false if surface was already marked as expired.
135  virtual bool expire();
136 
137  //- Update the surface as required.
138  // Do nothing (and return false) if no update was needed
139  virtual bool update();
140 
141 
142  //- Points of surface
143  virtual const pointField& points() const
144  {
145  return cuttingPlane::points();
146  }
147 
148  //- Faces of surface
149  virtual const faceList& faces() const
150  {
151  return cuttingPlane::faces();
152  }
153 
154  //- For every face original cell in mesh
155  const labelList& meshCells() const
156  {
157  return cuttingPlane::cutCells();
158  }
159 
160  //- Sample field on surface
161  virtual tmp<scalarField> sample
162  (
163  const volScalarField&
164  ) const;
165 
166 
167  //- Sample field on surface
169  (
170  const volVectorField&
171  ) const;
172 
173  //- Sample field on surface
175  (
177  ) const;
178 
179  //- Sample field on surface
181  (
182  const volSymmTensorField&
183  ) const;
184 
185  //- Sample field on surface
186  virtual tmp<tensorField> sample
187  (
188  const volTensorField&
189  ) const;
190 
191 
192  //- Interpolate field on surface
194  (
195  const interpolation<scalar>&
196  ) const;
197 
198 
199  //- Interpolate field on surface
201  (
202  const interpolation<vector>&
203  ) const;
204 
205  //- Interpolate field on surface
207  (
209  ) const;
210 
211  //- Interpolate field on surface
213  (
215  ) const;
216 
217  //- Interpolate field on surface
219  (
220  const interpolation<tensor>&
221  ) const;
222 
223  //- Write
224  virtual void print(Ostream&) const;
225 };
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace sampledSurfaces
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #ifdef NoRepository
236  #include "sampledPlaneTemplates.C"
237 #endif
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
dictionary dict
virtual const faceList & faces() const
Faces of surface.
virtual const faceList & faces() const
Faces of surface.
Definition: sampledPlane.H:168
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
TypeName("plane")
Runtime type information.
const word & name() const
Name of surface.
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPlane.C:88
bool interpolate() const
Interpolation requested for surface.
Generic GeometricField class.
virtual const pointField & points() const
Points of surface.
Definition: sampledPlane.H:162
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPlane.C:94
virtual ~plane()
Destructor.
Definition: sampledPlane.C:82
plane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledPlane.C:44
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledPlane.C:144
A class for handling words, derived from string.
Definition: word.H:59
virtual const pointField & points() const
Points of surface.
virtual bool update()
Update the surface as required.
Definition: sampledPlane.C:109
virtual void print(Ostream &) const
Write.
Definition: sampledPlane.C:232
const polyMesh & mesh() const
Access to the underlying mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const labelList & meshCells() const
For every face original cell in mesh.
Definition: sampledPlane.H:174
Abstract base class for interpolation.
const labelList & cutCells() const
Return List of cells cut by the plane.
Definition: cuttingPlane.H:160
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.