1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "sampledPlane.H"
27 #include "dictionary.H"
29 #include "volFields.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 defineTypeNameAndDebug(sampledPlane, 0);
38 addNamedToRunTimeSelectionTable(sampledSurface, sampledPlane, word, plane);
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 Foam::sampledPlane::sampledPlane
47 const plane& planeDesc,
51 sampledSurface(name, mesh),
52 cuttingPlane(planeDesc),
56 if (debug && zoneName_.size())
58 if (mesh.cellZones().findZoneID(zoneName_) < 0)
60 Info<< "cellZone \"" << zoneName_
61 << "\" not found - using entire mesh" << endl;
67 Foam::sampledPlane::sampledPlane
71 const dictionary& dict
74 sampledSurface(name, mesh, dict),
75 cuttingPlane(plane(dict.lookup("basePoint"), dict.lookup("normalVector"))),
76 zoneName_(word::null),
79 // make plane relative to the coordinateSystem (Cartesian)
80 // allow lookup from global coordinate systems
81 if (dict.found("coordinateSystem"))
83 coordinateSystem cs(dict, mesh);
85 point base = cs.globalPosition(planeDesc().refPoint());
86 vector norm = cs.globalVector(planeDesc().normal());
88 // assign the plane description
89 static_cast<plane&>(*this) = plane(base, norm);
92 dict.readIfPresent("zone", zoneName_);
94 if (debug && zoneName_.size())
96 if (mesh.cellZones().findZoneID(zoneName_) < 0)
98 Info<< "cellZone \"" << zoneName_
99 << "\" not found - using entire mesh" << endl;
106 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
108 Foam::sampledPlane::~sampledPlane()
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 bool Foam::sampledPlane::needsUpdate() const
120 bool Foam::sampledPlane::expire()
122 // already marked as expired
128 sampledSurface::clearGeom();
135 bool Foam::sampledPlane::update()
142 sampledSurface::clearGeom();
145 if (zoneName_.size())
147 zoneId = mesh().cellZones().findZoneID(zoneName_);
156 reCut(mesh(), mesh().cellZones()[zoneId]);
165 needsUpdate_ = false;
170 Foam::tmp<Foam::scalarField>
171 Foam::sampledPlane::sample
173 const volScalarField& vField
176 return sampleField(vField);
180 Foam::tmp<Foam::vectorField>
181 Foam::sampledPlane::sample
183 const volVectorField& vField
186 return sampleField(vField);
190 Foam::tmp<Foam::sphericalTensorField>
191 Foam::sampledPlane::sample
193 const volSphericalTensorField& vField
196 return sampleField(vField);
200 Foam::tmp<Foam::symmTensorField>
201 Foam::sampledPlane::sample
203 const volSymmTensorField& vField
206 return sampleField(vField);
210 Foam::tmp<Foam::tensorField>
211 Foam::sampledPlane::sample
213 const volTensorField& vField
216 return sampleField(vField);
220 Foam::tmp<Foam::scalarField>
221 Foam::sampledPlane::interpolate
223 const interpolation<scalar>& interpolator
226 return interpolateField(interpolator);
230 Foam::tmp<Foam::vectorField>
231 Foam::sampledPlane::interpolate
233 const interpolation<vector>& interpolator
236 return interpolateField(interpolator);
239 Foam::tmp<Foam::sphericalTensorField>
240 Foam::sampledPlane::interpolate
242 const interpolation<sphericalTensor>& interpolator
245 return interpolateField(interpolator);
249 Foam::tmp<Foam::symmTensorField>
250 Foam::sampledPlane::interpolate
252 const interpolation<symmTensor>& interpolator
255 return interpolateField(interpolator);
259 Foam::tmp<Foam::tensorField>
260 Foam::sampledPlane::interpolate
262 const interpolation<tensor>& interpolator
265 return interpolateField(interpolator);
269 void Foam::sampledPlane::print(Ostream& os) const
271 os << "sampledPlane: " << name() << " :"
272 << " base:" << refPoint()
273 << " normal:" << normal()
274 << " faces:" << faces().size()
275 << " points:" << points().size();
279 // ************************************************************************* //