1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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,
48 const keyType& zoneKey
51 sampledSurface(name, mesh),
52 cuttingPlane(planeDesc),
56 if (debug && zoneKey_.size() && mesh.cellZones().findIndex(zoneKey_) < 0)
58 Info<< "cellZone " << zoneKey_
59 << " not found - using entire mesh" << endl;
64 Foam::sampledPlane::sampledPlane
68 const dictionary& dict
71 sampledSurface(name, mesh, dict),
72 cuttingPlane(plane(dict.lookup("basePoint"), dict.lookup("normalVector"))),
73 zoneKey_(keyType::null),
76 // make plane relative to the coordinateSystem (Cartesian)
77 // allow lookup from global coordinate systems
78 if (dict.found("coordinateSystem"))
80 coordinateSystem cs(dict, mesh);
82 point base = cs.globalPosition(planeDesc().refPoint());
83 vector norm = cs.globalVector(planeDesc().normal());
85 // assign the plane description
86 static_cast<plane&>(*this) = plane(base, norm);
89 dict.readIfPresent("zone", zoneKey_);
91 if (debug && zoneKey_.size() && mesh.cellZones().findIndex(zoneKey_) < 0)
93 Info<< "cellZone " << zoneKey_
94 << " not found - using entire mesh" << endl;
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
101 Foam::sampledPlane::~sampledPlane()
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 bool Foam::sampledPlane::needsUpdate() const
113 bool Foam::sampledPlane::expire()
115 // already marked as expired
121 sampledSurface::clearGeom();
128 bool Foam::sampledPlane::update()
135 sampledSurface::clearGeom();
137 labelList selectedCells = mesh().cellZones().findMatching(zoneKey_).used();
139 if (selectedCells.empty())
141 reCut(mesh(), true); // always triangulate. Note:Make option?
145 reCut(mesh(), true, selectedCells);
154 needsUpdate_ = false;
159 Foam::tmp<Foam::scalarField> Foam::sampledPlane::sample
161 const volScalarField& vField
164 return sampleField(vField);
168 Foam::tmp<Foam::vectorField> Foam::sampledPlane::sample
170 const volVectorField& vField
173 return sampleField(vField);
177 Foam::tmp<Foam::sphericalTensorField> Foam::sampledPlane::sample
179 const volSphericalTensorField& vField
182 return sampleField(vField);
186 Foam::tmp<Foam::symmTensorField> Foam::sampledPlane::sample
188 const volSymmTensorField& vField
191 return sampleField(vField);
195 Foam::tmp<Foam::tensorField> Foam::sampledPlane::sample
197 const volTensorField& vField
200 return sampleField(vField);
204 Foam::tmp<Foam::scalarField> Foam::sampledPlane::interpolate
206 const interpolation<scalar>& interpolator
209 return interpolateField(interpolator);
213 Foam::tmp<Foam::vectorField> Foam::sampledPlane::interpolate
215 const interpolation<vector>& interpolator
218 return interpolateField(interpolator);
221 Foam::tmp<Foam::sphericalTensorField> Foam::sampledPlane::interpolate
223 const interpolation<sphericalTensor>& interpolator
226 return interpolateField(interpolator);
230 Foam::tmp<Foam::symmTensorField> Foam::sampledPlane::interpolate
232 const interpolation<symmTensor>& interpolator
235 return interpolateField(interpolator);
239 Foam::tmp<Foam::tensorField> Foam::sampledPlane::interpolate
241 const interpolation<tensor>& interpolator
244 return interpolateField(interpolator);
248 void Foam::sampledPlane::print(Ostream& os) const
250 os << "sampledPlane: " << name() << " :"
251 << " base:" << refPoint()
252 << " normal:" << normal()
253 << " faces:" << faces().size()
254 << " points:" << points().size();
258 // ************************************************************************* //