Loading...
Searching...
No Matches
sampledPatch.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "sampledPatch.H"
30#include "dictionary.H"
31#include "polyMesh.H"
32#include "polyPatch.H"
33#include "volFields.H"
34#include "surfaceFields.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
45}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const word& name,
53 const polyMesh& mesh,
55 const bool triangulate
56)
57:
59 selectionNames_(patchNames),
60 triangulate_(triangulate),
61 needsUpdate_(true)
62{}
63
64
66(
67 const word& name,
68 const polyMesh& mesh,
69 const dictionary& dict
70)
71:
73 selectionNames_(dict.get<wordRes>("patches")),
74 triangulate_(dict.getOrDefault("triangulate", false)),
75 needsUpdate_(true)
76{}
77
78
79// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80
82{
83 if (patchIDs_.empty())
84 {
85 labelList selected
86 (
87 mesh().boundaryMesh().patchSet(selectionNames_).sortedToc()
88 );
89
91 for (const label patchi : selected)
92 {
93 const polyPatch& pp = mesh().boundaryMesh()[patchi];
94
96 {
97 bad.append(patchi);
98 }
99 }
100
101 if (bad.size())
102 {
103 label nGood = (selected.size() - bad.size());
104
105 auto& os = nGood > 0 ? WarningInFunction : FatalErrorInFunction;
106
107 os << "Cannot sample an empty patch" << nl;
108
109 for (const label patchi : bad)
110 {
111 os << " "
112 << mesh().boundaryMesh()[patchi].name() << nl;
113 }
114
115 if (nGood)
116 {
117 os << "No non-empty patches selected" << endl
118 << exit(FatalError);
119 }
120 else
121 {
122 os << "Selected " << nGood << " non-empty patches" << nl;
123 }
124
125 patchIDs_.resize(nGood);
126 nGood = 0;
127 for (const label patchi : selected)
128 {
129 if (!bad.found(patchi))
130 {
131 patchIDs_[nGood] = patchi;
132 ++nGood;
133 }
134 }
135 }
136 else
137 {
138 patchIDs_ = std::move(selected);
140 }
141
142 return patchIDs_;
143}
144
147{
148 return needsUpdate_;
149}
150
151
153{
154 // already marked as expired
155 if (needsUpdate_)
156 {
157 return false;
158 }
159
161 Mesh::clear();
162
163 patchIDs_.clear();
164 patchStart_.clear();
165
166 patchIndex_.clear();
167 patchFaceLabels_.clear();
168
169 needsUpdate_ = true;
170 return true;
171}
172
173
175{
176 if (!needsUpdate_)
177 {
178 return false;
179 }
180
181 // Total number of faces selected
182 label numFaces = 0;
183 for (const label patchi : patchIDs())
184 {
185 const polyPatch& pp = mesh().boundaryMesh()[patchi];
186 numFaces += pp.size();
187 }
188
189 patchStart_.resize(patchIDs().size());
190
191 // The originating patch and local face in the patch.
192 patchIndex_.resize(numFaces);
193 patchFaceLabels_.resize(numFaces);
194
195 IndirectList<face> selectedFaces(mesh().faces(), labelList());
196 labelList& meshFaceIds = selectedFaces.addressing();
197 meshFaceIds.resize(numFaces);
198
199 numFaces = 0;
200
201 forAll(patchIDs(), idx)
202 {
203 const label patchi = patchIDs()[idx];
204 const polyPatch& pp = mesh().boundaryMesh()[patchi];
205 const label len = pp.size();
206
207 patchStart_[idx] = numFaces;
208
209 SubList<label>(patchIndex_, len, numFaces) = idx;
210
211 SubList<label>(patchFaceLabels_, len, numFaces) = identity(len);
212
213 SubList<label>(meshFaceIds, len, numFaces) = identity(len, pp.start());
214
215 numFaces += len;
216 }
217
218
219 uindirectPrimitivePatch allPatches(selectedFaces, mesh().points());
220
221 this->storedPoints() = allPatches.localPoints();
222 this->storedFaces() = allPatches.localFaces();
223
224
225 // triangulate uses remapFaces()
226 // - this is somewhat less efficient since it recopies the faces
227 // that we just created, but we probably don't want to do this
228 // too often anyhow.
229 if (triangulate_)
230 {
231 Mesh::triangulate();
232 }
233
234 if (debug)
235 {
236 print(Pout, debug);
237 Pout<< endl;
238 }
239
240 needsUpdate_ = false;
241 return true;
242}
243
244
245// remap action on triangulation
246void Foam::sampledPatch::remapFaces(const labelUList& faceMap)
247{
248 if (!faceMap.empty())
249 {
250 Mesh::remapFaces(faceMap);
251 patchFaceLabels_ = labelList
252 (
253 labelUIndList(patchFaceLabels_, faceMap)
254 );
255 patchIndex_ = labelList
256 (
257 labelUIndList(patchIndex_, faceMap)
258 );
259
260 // Update patchStart
261 if (patchIndex_.size())
262 {
263 patchStart_[patchIndex_[0]] = 0;
264 for (label i = 1; i < patchIndex_.size(); ++i)
265 {
266 if (patchIndex_[i] != patchIndex_[i-1])
267 {
268 patchStart_[patchIndex_[i]] = i;
270 }
271 }
272 }
273}
274
275
276Foam::tmp<Foam::scalarField> Foam::sampledPatch::sample
277(
279) const
280{
281 return sampleOnFaces(sampler);
282}
283
284
286(
288) const
289{
290 return sampleOnFaces(sampler);
291}
292
293
295(
297) const
298{
299 return sampleOnFaces(sampler);
300}
301
302
304(
306) const
307{
308 return sampleOnFaces(sampler);
309}
310
311
313(
315) const
316{
317 return sampleOnFaces(sampler);
318}
319
322{
323 return true;
324}
325
326
328(
329 const surfaceScalarField& sField
330) const
331{
332 return sampleOnFaces(sField);
333}
334
335
337(
338 const surfaceVectorField& sField
339) const
340{
341 return sampleOnFaces(sField);
342}
343
344
346(
347 const surfaceSphericalTensorField& sField
348) const
349{
350 return sampleOnFaces(sField);
351}
352
353
355(
356 const surfaceSymmTensorField& sField
357) const
358{
359 return sampleOnFaces(sField);
360}
361
362
363Foam::tmp<Foam::tensorField> Foam::sampledPatch::sample
364(
365 const surfaceTensorField& sField
366) const
367{
368 return sampleOnFaces(sField);
369}
370
371
373(
374 const interpolation<scalar>& interpolator
375) const
376{
377 return sampleOnPoints(interpolator);
378}
379
380
382(
383 const interpolation<vector>& interpolator
384) const
385{
386 return sampleOnPoints(interpolator);
387}
388
389
391(
393) const
394{
395 return sampleOnPoints(interpolator);
396}
397
398
400(
401 const interpolation<symmTensor>& interpolator
402) const
403{
404 return sampleOnPoints(interpolator);
405}
406
407
409(
410 const interpolation<tensor>& interpolator
411) const
412{
413 return sampleOnPoints(interpolator);
414}
415
416
417void Foam::sampledPatch::print(Ostream& os, int level) const
418{
419 os << "sampledPatch: " << name() << " :"
420 << " patches:" << flatOutput(selectionNames_);
421
422 if (level)
423 {
424 os << " faces:" << faces().size()
425 << " points:" << points().size();
426 }
427}
428
429
430// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
label size() const noexcept
The number of elements in the list.
A List with indirect addressing.
const Addr & addressing() const noexcept
The list addressing.
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
virtual label triangulate()
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
bool get(const label i) const
Definition UList.H:868
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for volume field interpolation.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A sampledSurface on patches. Non-triangulated by default.
virtual void print(Ostream &os, int level=0) const
Print information.
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample boundary of volume field onto surface faces.
const labelList & patchIDs() const
The patches selected.
virtual const faceList & faces() const
Faces of surface.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool withSurfaceFields() const
Can it sample surface-fields?
virtual bool update()
Update the surface as required.
const wordRes & patchNames() const
The selection (word/regex) of patches.
sampledPatch(const word &name, const polyMesh &mesh, const UList< wordRe > &patchNames, const bool triangulate=false)
Construct from components.
An abstract class for surfaces with sampling.
sampledSurface(const word &name, std::nullptr_t)
Construct null.
const word & name() const noexcept
Name of surface.
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
bool interpolate() const noexcept
Same as isPointData().
A class for managing temporary objects.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
Ostream & print(Ostream &os, UIntType value, char off='0', char on='1')
Print 0/1 bits in the (unsigned) integral type.
Definition BitOps.H:320
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< label > labelList
A List of labels.
Definition List.H:62
UIndirectList< label > labelUIndList
UIndirectList of labels.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
GeometricField< tensor, fvsPatchField, surfaceMesh > surfaceTensorField
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
GeometricField< sphericalTensor, fvsPatchField, surfaceMesh > surfaceSphericalTensorField
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
GeometricField< symmTensor, fvsPatchField, surfaceMesh > surfaceSymmTensorField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.