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 "sampledPatch.H"
28 #include "dictionary.H"
30 #include "polyPatch.H"
31 #include "volFields.H"
33 #include "addToRunTimeSelectionTable.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(sampledPatch, 0);
40 addNamedToRunTimeSelectionTable(sampledSurface, sampledPatch, word, patch);
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 Foam::sampledPatch::sampledPatch
50 const word& patchName,
51 const bool triangulate
54 sampledSurface(name, mesh),
55 patchName_(patchName),
56 triangulate_(triangulate),
62 Foam::sampledPatch::sampledPatch
66 const dictionary& dict
69 sampledSurface(name, mesh, dict),
70 patchName_(dict.lookup("patchName")),
71 triangulate_(dict.lookupOrDefault("triangulate", false)),
77 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 Foam::sampledPatch::~sampledPatch()
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 bool Foam::sampledPatch::needsUpdate() const
91 bool Foam::sampledPatch::expire()
93 // already marked as expired
99 sampledSurface::clearGeom();
100 MeshStorage::clear();
101 patchFaceLabels_.clear();
108 bool Foam::sampledPatch::update()
115 label patchI = mesh().boundaryMesh().findPatchID(patchName_);
119 const polyPatch& p = mesh().boundaryMesh()[patchI];
120 this->storedPoints() = p.localPoints();
121 this->storedFaces() = p.localFaces();
124 patchFaceLabels_.setSize(faces().size());
125 forAll(patchFaceLabels_, i)
127 patchFaceLabels_[i] = i;
130 // triangulate uses remapFaces()
131 // - this is somewhat less efficient since it recopies the faces
132 // that we just created, but we probably don't want to do this
136 MeshStorage::triangulate();
146 needsUpdate_ = false;
151 // remap action on triangulation
152 void Foam::sampledPatch::remapFaces
154 const UList<label>& faceMap
157 // recalculate the cells cut
158 if (&faceMap && faceMap.size())
160 MeshStorage::remapFaces(faceMap);
166 Foam::tmp<Foam::scalarField>
167 Foam::sampledPatch::sample
169 const volScalarField& vField
172 return sampleField(vField);
176 Foam::tmp<Foam::vectorField>
177 Foam::sampledPatch::sample
179 const volVectorField& vField
182 return sampleField(vField);
185 Foam::tmp<Foam::sphericalTensorField>
186 Foam::sampledPatch::sample
188 const volSphericalTensorField& vField
191 return sampleField(vField);
195 Foam::tmp<Foam::symmTensorField>
196 Foam::sampledPatch::sample
198 const volSymmTensorField& vField
201 return sampleField(vField);
205 Foam::tmp<Foam::tensorField>
206 Foam::sampledPatch::sample
208 const volTensorField& vField
211 return sampleField(vField);
215 Foam::tmp<Foam::scalarField>
216 Foam::sampledPatch::interpolate
218 const interpolation<scalar>& interpolator
221 return interpolateField(interpolator);
225 Foam::tmp<Foam::vectorField>
226 Foam::sampledPatch::interpolate
228 const interpolation<vector>& interpolator
231 return interpolateField(interpolator);
234 Foam::tmp<Foam::sphericalTensorField>
235 Foam::sampledPatch::interpolate
237 const interpolation<sphericalTensor>& interpolator
240 return interpolateField(interpolator);
244 Foam::tmp<Foam::symmTensorField>
245 Foam::sampledPatch::interpolate
247 const interpolation<symmTensor>& interpolator
250 return interpolateField(interpolator);
254 Foam::tmp<Foam::tensorField>
255 Foam::sampledPatch::interpolate
257 const interpolation<tensor>& interpolator
260 return interpolateField(interpolator);
264 void Foam::sampledPatch::print(Ostream& os) const
266 os << "sampledPatch: " << name() << " :"
267 << " patch:" << patchName()
268 << " faces:" << faces().size()
269 << " points:" << points().size();
273 // ************************************************************************* //