sampledPlane.C
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-2020 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 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledPlane.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace sampledSurfaces
34 {
35  defineTypeNameAndDebug(plane, 0);
36  addToRunTimeSelectionTable(sampledSurface, plane, word);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const word& name,
46  const polyMesh& mesh,
47  const dictionary& dict
48 )
49 :
50  sampledSurface(name, mesh, dict),
52  zoneKey_(dict.lookupOrDefault("zone", wordRe::null)),
53  triangulate_(dict.lookupOrDefault("triangulate", true)),
54  needsUpdate_(true)
55 {
56  // Make plane relative to the coordinateSystem (Cartesian)
57  // allow lookup from global coordinate systems
58  if (dict.found("coordinateSystem"))
59  {
61  (
62  coordinateSystem::New(mesh, dict.subDict("coordinateSystem"))
63  );
64 
65  point base = cs->globalPosition(planeDesc().refPoint());
66  vector norm = cs->globalVector(planeDesc().normal());
67 
68  // Assign the plane description
69  static_cast<Foam::plane&>(*this) = Foam::plane(base, norm);
70  }
71 
72  if (debug && zoneKey_.size() && mesh.cellZones().findIndex(zoneKey_) < 0)
73  {
74  Info<< "cellZone " << zoneKey_
75  << " not found - using entire mesh" << endl;
76  }
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {
90  return needsUpdate_;
91 }
92 
93 
95 {
96  // Already marked as expired
97  if (needsUpdate_)
98  {
99  return false;
100  }
101 
103 
104  needsUpdate_ = true;
105  return true;
106 }
107 
108 
110 {
111  if (!needsUpdate_)
112  {
113  return false;
114  }
115 
117 
118  const labelList selectedCells
119  (
120  mesh().cellZones().findMatching(zoneKey_).used()
121  );
122 
123  if (selectedCells.empty())
124  {
125  reCut(mesh(), triangulate_);
126  }
127  else
128  {
129  reCut(mesh(), triangulate_, selectedCells);
130  }
131 
132  if (debug)
133  {
134  print(Pout);
135  Pout<< endl;
136  }
137 
138  needsUpdate_ = false;
139  return true;
140 }
141 
142 
144 (
145  const volScalarField& vField
146 ) const
147 {
148  return sampleField(vField);
149 }
150 
151 
153 (
154  const volVectorField& vField
155 ) const
156 {
157  return sampleField(vField);
158 }
159 
160 
162 (
163  const volSphericalTensorField& vField
164 ) const
165 {
166  return sampleField(vField);
167 }
168 
169 
171 (
172  const volSymmTensorField& vField
173 ) const
174 {
175  return sampleField(vField);
176 }
177 
178 
180 (
181  const volTensorField& vField
182 ) const
183 {
184  return sampleField(vField);
185 }
186 
187 
189 (
190  const interpolation<scalar>& interpolator
191 ) const
192 {
193  return interpolateField(interpolator);
194 }
195 
196 
198 (
199  const interpolation<vector>& interpolator
200 ) const
201 {
202  return interpolateField(interpolator);
203 }
204 
206 (
207  const interpolation<sphericalTensor>& interpolator
208 ) const
209 {
210  return interpolateField(interpolator);
211 }
212 
213 
215 (
216  const interpolation<symmTensor>& interpolator
217 ) const
218 {
219  return interpolateField(interpolator);
220 }
221 
222 
224 (
225  const interpolation<tensor>& interpolator
226 ) const
227 {
228  return interpolateField(interpolator);
229 }
230 
231 
233 {
234  os << "plane: " << name() << " :"
235  << " base:" << refPoint()
236  << " normal:" << normal()
237  << " triangulate:" << triangulate_
238  << " faces:" << faces().size()
239  << " points:" << points().size();
240 }
241 
242 
243 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
An abstract class for surfaces with sampling.
defineTypeNameAndDebug(distanceSurface, 0)
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPlane.C:88
bool interpolate() const
Interpolation requested for surface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
fvMesh & mesh
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPlane.C:94
Macros for easy insertion into run-time selection tables.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
A sampledSurface defined by a plane using the iso-surface algorithm to &#39;cut&#39; the mesh.
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 void clearGeom() const
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
virtual bool update()
Update the surface as required.
Definition: sampledPlane.C:109
virtual void print(Ostream &) const
Write.
Definition: sampledPlane.C:232
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledPlane.C:144
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static autoPtr< coordinateSystem > New(const objectRegistry &obr, const dictionary &dict)
Select constructed from dictionary and objectRegistry.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const wordRe &) const
Return zone index for the first match, return -1 if not found.
Definition: MeshZones.C:306
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A class for managing temporary objects.
Definition: PtrList.H:53
static const wordRe null
An empty wordRe.
Definition: wordRe.H:88
Namespace for OpenFOAM.