sampledPlane.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledPlane.H"
27 #include "dictionary.H"
28 #include "polyMesh.H"
29 #include "volFields.H"
30 
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(sampledPlane, 0);
38  addNamedToRunTimeSelectionTable(sampledSurface, sampledPlane, word, plane);
39 }
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const word& name,
46  const polyMesh& mesh,
47  const plane& planeDesc,
48  const keyType& zoneKey,
49  const bool triangulate
50 )
51 :
53  cuttingPlane(planeDesc),
54  zoneKey_(zoneKey),
55  triangulate_(triangulate),
56  needsUpdate_(true)
57 {
58  if (debug && zoneKey_.size() && mesh.cellZones().findIndex(zoneKey_) < 0)
59  {
60  Info<< "cellZone " << zoneKey_
61  << " not found - using entire mesh" << endl;
62  }
63 }
64 
65 
67 (
68  const word& name,
69  const polyMesh& mesh,
70  const dictionary& dict
71 )
72 :
74  cuttingPlane(plane(dict.lookup("basePoint"), dict.lookup("normalVector"))),
75  zoneKey_(keyType::null),
76  triangulate_(dict.lookupOrDefault("triangulate", true)),
77  needsUpdate_(true)
78 {
79  // make plane relative to the coordinateSystem (Cartesian)
80  // allow lookup from global coordinate systems
81  if (dict.found("coordinateSystem"))
82  {
83  coordinateSystem cs(mesh, dict.subDict("coordinateSystem"));
84 
85  point base = cs.globalPosition(planeDesc().refPoint());
86  vector norm = cs.globalVector(planeDesc().normal());
87 
88  // assign the plane description
89  static_cast<plane&>(*this) = plane(base, norm);
90  }
91 
92  dict.readIfPresent("zone", zoneKey_);
93 
94  if (debug && zoneKey_.size() && mesh.cellZones().findIndex(zoneKey_) < 0)
95  {
96  Info<< "cellZone " << zoneKey_
97  << " not found - using entire mesh" << endl;
98  }
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
103 
105 {}
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
111 {
112  return needsUpdate_;
113 }
114 
115 
117 {
118  // already marked as expired
119  if (needsUpdate_)
120  {
121  return false;
122  }
123 
125 
126  needsUpdate_ = true;
127  return true;
128 }
129 
130 
132 {
133  if (!needsUpdate_)
134  {
135  return false;
136  }
137 
139 
140  labelList selectedCells = mesh().cellZones().findMatching(zoneKey_).used();
141 
142  if (selectedCells.empty())
143  {
144  reCut(mesh(), triangulate_);
145  }
146  else
147  {
148  reCut(mesh(), triangulate_, selectedCells);
149  }
150 
151  if (debug)
152  {
153  print(Pout);
154  Pout<< endl;
155  }
156 
157  needsUpdate_ = false;
158  return true;
159 }
160 
161 
163 (
164  const volScalarField& vField
165 ) const
166 {
167  return sampleField(vField);
168 }
169 
170 
172 (
173  const volVectorField& vField
174 ) const
175 {
176  return sampleField(vField);
177 }
178 
179 
181 (
182  const volSphericalTensorField& vField
183 ) const
184 {
185  return sampleField(vField);
186 }
187 
188 
190 (
191  const volSymmTensorField& vField
192 ) const
193 {
194  return sampleField(vField);
195 }
196 
197 
199 (
200  const volTensorField& vField
201 ) const
202 {
203  return sampleField(vField);
204 }
205 
206 
208 (
209  const interpolation<scalar>& interpolator
210 ) const
211 {
212  return interpolateField(interpolator);
213 }
214 
215 
217 (
218  const interpolation<vector>& interpolator
219 ) const
220 {
221  return interpolateField(interpolator);
222 }
223 
225 (
226  const interpolation<sphericalTensor>& interpolator
227 ) const
228 {
229  return interpolateField(interpolator);
230 }
231 
232 
234 (
235  const interpolation<symmTensor>& interpolator
236 ) const
237 {
238  return interpolateField(interpolator);
239 }
240 
241 
243 (
244  const interpolation<tensor>& interpolator
245 ) const
246 {
247  return interpolateField(interpolator);
248 }
249 
250 
252 {
253  os << "sampledPlane: " << name() << " :"
254  << " base:" << refPoint()
255  << " normal:" << normal()
256  << " triangulate:" << triangulate_
257  << " faces:" << faces().size()
258  << " points:" << points().size();
259 }
260 
261 
262 // ************************************************************************* //
Foam::coordinateSystem::globalPosition
point globalPosition(const point &local) const
Convert from position in local coordinate system to global.
Definition: coordinateSystem.H:320
volFields.H
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch, cyclicPeriodicAMI)
Foam::sampledPlane::update
virtual bool update()
Update the surface as required.
Definition: sampledPlane.C:131
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
print
void print(const char *msg, Ostream &os, const PtrList< GeoField > &flds)
Definition: foamToVTK.C:164
Foam::sampledPlane::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPlane.C:116
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::cuttingPlane
Constructs plane through mesh.
Definition: cuttingPlane.H:60
Foam::sampledPlane::sampledPlane
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
Foam::sampledPlane::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledPlane.C:163
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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::keyType
A class for handling keywords in dictionaries.
Definition: keyType.H:56
sampledPlane.H
Foam::sampledPlane::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPlane.C:110
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::Info
messageStream Info
Foam::coordinateSystem::globalVector
vector globalVector(const vector &local) const
Convert from vector components in local coordinate system to.
Definition: coordinateSystem.H:334
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::interpolation< scalar >
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::ZoneMesh::findMatching
PackedBoolList findMatching(const keyType &) const
Mark cells that match the zone specification.
Definition: ZoneMesh.C:377
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::sampledPlane::print
virtual void print(Ostream &) const
Write.
Definition: sampledPlane.C:251
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::Vector< scalar >
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
dictionary.H
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Definition: sampledSurface.C:42
Foam::sampledPlane::~sampledPlane
virtual ~sampledPlane()
Destructor.
Definition: sampledPlane.C:104
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::sampledSurface::interpolate
bool interpolate() const
Interpolation requested for surface.
Definition: sampledSurface.H:256
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::PackedBoolList::used
Xfer< labelList > used() const
Return indices of the used (true) elements as a list of labels.
Definition: PackedBoolList.C:264
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
normal
A normal distribution model.
Foam::coordinateSystem
Base class for other coordinate system specifications.
Definition: coordinateSystem.H:85