1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "sampledPlane.H"
28 #include "dictionary.H"
30 #include "volFields.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(sampledPlane, 0);
39 addNamedToRunTimeSelectionTable(sampledSurface, sampledPlane, word, plane);
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 Foam::sampledPlane::sampledPlane
48 const plane& planeDesc,
52 sampledSurface(name, mesh),
53 cuttingPlane(planeDesc),
57 if (debug && zoneName_.size())
59 if (mesh.cellZones().findZoneID(zoneName_) < 0)
61 Info<< "cellZone \"" << zoneName_
62 << "\" not found - using entire mesh" << endl;
68 Foam::sampledPlane::sampledPlane
72 const dictionary& dict
75 sampledSurface(name, mesh, dict),
76 cuttingPlane(plane(dict.lookup("basePoint"), dict.lookup("normalVector"))),
77 zoneName_(word::null),
80 // make plane relative to the coordinateSystem (Cartesian)
81 // allow lookup from global coordinate systems
82 if (dict.found("coordinateSystem"))
84 coordinateSystem cs(dict, mesh);
86 point base = cs.globalPosition(planeDesc().refPoint());
87 vector norm = cs.globalVector(planeDesc().normal());
89 // assign the plane description
90 static_cast<plane&>(*this) = plane(base, norm);
93 dict.readIfPresent("zone", zoneName_);
95 if (debug && zoneName_.size())
97 if (mesh.cellZones().findZoneID(zoneName_) < 0)
99 Info<< "cellZone \"" << zoneName_
100 << "\" not found - using entire mesh" << endl;
107 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
109 Foam::sampledPlane::~sampledPlane()
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 bool Foam::sampledPlane::needsUpdate() const
121 bool Foam::sampledPlane::expire()
123 // already marked as expired
129 sampledSurface::clearGeom();
136 bool Foam::sampledPlane::update()
143 sampledSurface::clearGeom();
146 if (zoneName_.size())
148 zoneId = mesh().cellZones().findZoneID(zoneName_);
157 reCut(mesh(), mesh().cellZones()[zoneId]);
166 needsUpdate_ = false;
171 Foam::tmp<Foam::scalarField>
172 Foam::sampledPlane::sample
174 const volScalarField& vField
177 return sampleField(vField);
181 Foam::tmp<Foam::vectorField>
182 Foam::sampledPlane::sample
184 const volVectorField& vField
187 return sampleField(vField);
191 Foam::tmp<Foam::sphericalTensorField>
192 Foam::sampledPlane::sample
194 const volSphericalTensorField& vField
197 return sampleField(vField);
201 Foam::tmp<Foam::symmTensorField>
202 Foam::sampledPlane::sample
204 const volSymmTensorField& vField
207 return sampleField(vField);
211 Foam::tmp<Foam::tensorField>
212 Foam::sampledPlane::sample
214 const volTensorField& vField
217 return sampleField(vField);
221 Foam::tmp<Foam::scalarField>
222 Foam::sampledPlane::interpolate
224 const interpolation<scalar>& interpolator
227 return interpolateField(interpolator);
231 Foam::tmp<Foam::vectorField>
232 Foam::sampledPlane::interpolate
234 const interpolation<vector>& interpolator
237 return interpolateField(interpolator);
240 Foam::tmp<Foam::sphericalTensorField>
241 Foam::sampledPlane::interpolate
243 const interpolation<sphericalTensor>& interpolator
246 return interpolateField(interpolator);
250 Foam::tmp<Foam::symmTensorField>
251 Foam::sampledPlane::interpolate
253 const interpolation<symmTensor>& interpolator
256 return interpolateField(interpolator);
260 Foam::tmp<Foam::tensorField>
261 Foam::sampledPlane::interpolate
263 const interpolation<tensor>& interpolator
266 return interpolateField(interpolator);
270 void Foam::sampledPlane::print(Ostream& os) const
272 os << "sampledPlane: " << name() << " :"
273 << " base:" << refPoint()
274 << " normal:" << normal()
275 << " faces:" << faces().size()
276 << " points:" << points().size();
280 // ************************************************************************* //