sampledCuttingPlane.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-2015 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::sampledCuttingPlane
26 
27 Description
28  A sampledSurface defined by a plane
29 
30 SourceFiles
31  sampledCuttingPlane.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef sampledCuttingPlane_H
36 #define sampledCuttingPlane_H
37 
38 #include "sampledSurface.H"
39 #include "isoSurface.H"
40 //#include "isoSurfaceCell.H"
41 #include "plane.H"
42 #include "ZoneIDs.H"
43 #include "fvMeshSubset.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class sampledCuttingPlane Declaration
52 \*---------------------------------------------------------------------------*/
53 
55 :
56  public sampledSurface
57 {
58  // Private data
59 
60  //- Plane
61  const plane plane_;
62 
63  //- Optional bounding box to trim triangles against
64  const boundBox bounds_;
65 
66  //- Merge tolerance
67  const scalar mergeTol_;
68 
69  //- Whether to coarsen
70  const Switch regularise_;
71 
72  //- Whether to recalculate cell values as average of point values
73  const Switch average_;
74 
75  //- Zone name/index (if restricted to zones)
76  mutable cellZoneID zoneID_;
77 
78  //- For zones: patch to put exposed faces into
79  mutable word exposedPatchName_;
80 
81  //- Track if the surface needs an update
82  mutable bool needsUpdate_;
83 
84 
85  //- Optional subsetted mesh
87 
88  //- Distance to cell centres
90 
91  //- Distance to points
93 
94  //- Constructed iso surface
95  //autoPtr<isoSurfaceCell> isoSurfPtr_;
97 
98  //- Triangles converted to faceList
100 
101 
102  // Private Member Functions
103 
104  //- Create iso surface
105  void createGeometry();
106 
107  //- Sample field on faces
108  template<class Type>
110  (
112  ) const;
113 
114 
115  template<class Type>
118 
119 
120 public:
121 
122  //- Runtime type information
123  TypeName("sampledCuttingPlane");
124 
125 
126  // Constructors
127 
128  //- Construct from dictionary
130  (
131  const word& name,
132  const polyMesh& mesh,
133  const dictionary& dict
134  );
135 
136 
137  //- Destructor
138  virtual ~sampledCuttingPlane();
139 
140 
141  // Member Functions
142 
143  //- Does the surface need an update?
144  virtual bool needsUpdate() const;
145 
146  //- Mark the surface as needing an update.
147  // May also free up unneeded data.
148  // Return false if surface was already marked as expired.
149  virtual bool expire();
150 
151  //- Update the surface as required.
152  // Do nothing (and return false) if no update was needed
153  virtual bool update();
154 
155  //- Points of surface
156  virtual const pointField& points() const
157  {
158  return surface().points();
159  }
160 
161  //- Faces of surface
162  virtual const faceList& faces() const
163  {
164  if (facesPtr_.empty())
165  {
166  const triSurface& s = surface();
167 
168  facesPtr_.reset(new faceList(s.size()));
169 
170  forAll(s, i)
171  {
172  facesPtr_()[i] = s[i].triFaceFace();
173  }
174  }
175  return facesPtr_;
176  }
177 
178 
179  //const isoSurfaceCell& surface() const
180  const isoSurface& surface() const
181  {
182  return isoSurfPtr_();
183  }
184 
185  //- Sample field on surface
186  virtual tmp<scalarField> sample
187  (
188  const volScalarField&
189  ) const;
190 
191  //- Sample field on surface
192  virtual tmp<vectorField> sample
193  (
194  const volVectorField&
195  ) const;
196 
197  //- Sample field on surface
199  (
201  ) const;
202 
203  //- Sample field on surface
205  (
206  const volSymmTensorField&
207  ) const;
208 
209  //- Sample field on surface
210  virtual tmp<tensorField> sample
211  (
212  const volTensorField&
213  ) const;
214 
215 
216  //- Interpolate field on surface
218  (
219  const interpolation<scalar>&
220  ) const;
221 
222  //- Interpolate field on surface
224  (
225  const interpolation<vector>&
226  ) const;
227 
228  //- Interpolate field on surface
230  (
232  ) const;
233 
234  //- Interpolate field on surface
236  (
238  ) const;
239 
240  //- Interpolate field on surface
242  (
243  const interpolation<tensor>&
244  ) const;
245 
246  //- Write
247  virtual void print(Ostream&) const;
248 };
249 
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 } // End namespace Foam
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 #ifdef NoRepository
259 #endif
260 
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263 #endif
264 
265 // ************************************************************************* //
Foam::sampledCuttingPlane::points
virtual const pointField & points() const
Points of surface.
Definition: sampledCuttingPlane.H:155
Foam::sampledCuttingPlane::cellDistancePtr_
autoPtr< volScalarField > cellDistancePtr_
Distance to cell centres.
Definition: sampledCuttingPlane.H:88
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::sampledCuttingPlane::isoSurfPtr_
autoPtr< isoSurface > isoSurfPtr_
Constructed iso surface.
Definition: sampledCuttingPlane.H:95
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::sampledCuttingPlane::average_
const Switch average_
Whether to recalculate cell values as average of point values.
Definition: sampledCuttingPlane.H:72
Foam::sampledCuttingPlane::TypeName
TypeName("sampledCuttingPlane")
Runtime type information.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::sampledCuttingPlane
A sampledSurface defined by a plane.
Definition: sampledCuttingPlane.H:53
sampledCuttingPlaneTemplates.C
fvMeshSubset.H
ZoneIDs.H
Foam::sampledCuttingPlane::surface
const isoSurface & surface() const
Definition: sampledCuttingPlane.H:179
Foam::DynamicID< cellZoneMesh >
Foam::sampledCuttingPlane::regularise_
const Switch regularise_
Whether to coarsen.
Definition: sampledCuttingPlane.H:69
Foam::sampledCuttingPlane::print
virtual void print(Ostream &) const
Write.
Definition: sampledCuttingPlane.C:461
Foam::sampledCuttingPlane::createGeometry
void createGeometry()
Create iso surface.
Definition: sampledCuttingPlane.C:49
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::sampledCuttingPlane::plane_
const plane plane_
Plane.
Definition: sampledCuttingPlane.H:60
Foam::plane
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
Foam::sampledCuttingPlane::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledCuttingPlane.C:309
Foam::sampledCuttingPlane::bounds_
const boundBox bounds_
Optional bounding box to trim triangles against.
Definition: sampledCuttingPlane.H:63
Foam::isoSurface
A surface formed by the iso value. After "Regularised Marching Tetrahedra: improved iso-surface extra...
Definition: isoSurface.H:85
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
plane.H
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
sampledSurface.H
Foam::sampledCuttingPlane::mergeTol_
const scalar mergeTol_
Merge tolerance.
Definition: sampledCuttingPlane.H:66
Foam::sampledCuttingPlane::update
virtual bool update()
Update the surface as required.
Definition: sampledCuttingPlane.C:341
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::sampledCuttingPlane::~sampledCuttingPlane
virtual ~sampledCuttingPlane()
Destructor.
Definition: sampledCuttingPlane.C:303
Foam::interpolation< Type >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::sampledCuttingPlane::interpolateField
tmp< Field< Type > > interpolateField(const interpolation< Type > &) const
Foam::sampledCuttingPlane::subMeshPtr_
autoPtr< fvMeshSubset > subMeshPtr_
Optional subsetted mesh.
Definition: sampledCuttingPlane.H:85
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::sampledCuttingPlane::needsUpdate_
bool needsUpdate_
Track if the surface needs an update.
Definition: sampledCuttingPlane.H:81
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::sampledCuttingPlane::sampledCuttingPlane
sampledCuttingPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledCuttingPlane.C:258
Foam::sampledCuttingPlane::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledCuttingPlane.C:364
isoSurface.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::sampledCuttingPlane::zoneID_
cellZoneID zoneID_
Zone name/index (if restricted to zones)
Definition: sampledCuttingPlane.H:75
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sampledCuttingPlane::pointDistance_
scalarField pointDistance_
Distance to points.
Definition: sampledCuttingPlane.H:91
Foam::sampledSurface::name
const word & name() const
Name of surface.
Definition: sampledSurface.H:250
Foam::sampledCuttingPlane::faces
virtual const faceList & faces() const
Faces of surface.
Definition: sampledCuttingPlane.H:161
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:244
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::sampledCuttingPlane::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledCuttingPlane.C:315
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation requested for surface.
Definition: sampledSurface.H:256
Foam::sampledCuttingPlane::facesPtr_
autoPtr< faceList > facesPtr_
Triangles converted to faceList.
Definition: sampledCuttingPlane.H:98
Foam::sampledCuttingPlane::exposedPatchName_
word exposedPatchName_
For zones: patch to put exposed faces into.
Definition: sampledCuttingPlane.H:78
Foam::sampledCuttingPlane::sampleField
tmp< Field< Type > > sampleField(const GeometricField< Type, fvPatchField, volMesh > &vField) const
Sample field on faces.