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 "sampledPatch.H"
27 #include "dictionary.H"
29 #include "polyPatch.H"
30 #include "volFields.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(sampledPatch, 0);
39 addNamedToRunTimeSelectionTable(sampledSurface, sampledPatch, word, patch);
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 Foam::sampledPatch::sampledPatch
49 const word& patchName,
50 const bool triangulate
53 sampledSurface(name, mesh),
54 patchName_(patchName),
55 triangulate_(triangulate),
61 Foam::sampledPatch::sampledPatch
65 const dictionary& dict
68 sampledSurface(name, mesh, dict),
69 patchName_(dict.lookup("patchName")),
70 triangulate_(dict.lookupOrDefault("triangulate", false)),
76 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
78 Foam::sampledPatch::~sampledPatch()
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 bool Foam::sampledPatch::needsUpdate() const
90 bool Foam::sampledPatch::expire()
92 // already marked as expired
98 sampledSurface::clearGeom();
100 patchFaceLabels_.clear();
107 bool Foam::sampledPatch::update()
114 const label patchI = mesh().boundaryMesh().findPatchID(patchName_);
118 const polyPatch& p = mesh().boundaryMesh()[patchI];
119 this->storedPoints() = p.localPoints();
120 this->storedFaces() = p.localFaces();
123 patchFaceLabels_.setSize(faces().size());
124 forAll(patchFaceLabels_, i)
126 patchFaceLabels_[i] = i;
129 // triangulate uses remapFaces()
130 // - this is somewhat less efficient since it recopies the faces
131 // that we just created, but we probably don't want to do this
135 MeshStorage::triangulate();
145 needsUpdate_ = false;
150 // remap action on triangulation
151 void Foam::sampledPatch::remapFaces
153 const UList<label>& faceMap
156 // recalculate the cells cut
157 if (&faceMap && faceMap.size())
159 MeshStorage::remapFaces(faceMap);
160 patchFaceLabels_ = labelList
162 UIndirectList<label>(patchFaceLabels_, faceMap)
169 Foam::tmp<Foam::scalarField>
170 Foam::sampledPatch::sample
172 const volScalarField& vField
175 return sampleField(vField);
179 Foam::tmp<Foam::vectorField>
180 Foam::sampledPatch::sample
182 const volVectorField& vField
185 return sampleField(vField);
188 Foam::tmp<Foam::sphericalTensorField>
189 Foam::sampledPatch::sample
191 const volSphericalTensorField& vField
194 return sampleField(vField);
198 Foam::tmp<Foam::symmTensorField>
199 Foam::sampledPatch::sample
201 const volSymmTensorField& vField
204 return sampleField(vField);
208 Foam::tmp<Foam::tensorField>
209 Foam::sampledPatch::sample
211 const volTensorField& vField
214 return sampleField(vField);
218 Foam::tmp<Foam::scalarField>
219 Foam::sampledPatch::interpolate
221 const interpolation<scalar>& interpolator
224 return interpolateField(interpolator);
228 Foam::tmp<Foam::vectorField>
229 Foam::sampledPatch::interpolate
231 const interpolation<vector>& interpolator
234 return interpolateField(interpolator);
237 Foam::tmp<Foam::sphericalTensorField>
238 Foam::sampledPatch::interpolate
240 const interpolation<sphericalTensor>& interpolator
243 return interpolateField(interpolator);
247 Foam::tmp<Foam::symmTensorField>
248 Foam::sampledPatch::interpolate
250 const interpolation<symmTensor>& interpolator
253 return interpolateField(interpolator);
257 Foam::tmp<Foam::tensorField>
258 Foam::sampledPatch::interpolate
260 const interpolation<tensor>& interpolator
263 return interpolateField(interpolator);
267 void Foam::sampledPatch::print(Ostream& os) const
269 os << "sampledPatch: " << name() << " :"
270 << " patch:" << patchName()
271 << " faces:" << faces().size()
272 << " points:" << points().size();
276 // ************************************************************************* //