Loading...
Searching...
No Matches
polyPatch.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-2025 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 "polyPatch.H"
31#include "polyBoundaryMesh.H"
32#include "polyMesh.H"
33#include "primitiveMesh.H"
34#include "SubField.H"
35#include "entry.H"
36#include "dictionary.H"
37#include "pointPatchField.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44
46 (
47 debug::debugSwitch("disallowGenericPolyPatch", 0)
48 );
49
52
57
58// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
64
65
70}
71
72
76}
77
78
79// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80
82(
83 const word& name,
84 const label size,
85 const label start,
86 const label index,
87 const polyBoundaryMesh& bm,
88 const word& patchType
89)
90:
91 patchIdentifier(name, index),
93 (
94 SubList<face>(bm.mesh().faces(), size, start),
95 bm.mesh().points()
96 ),
97 start_(start),
98 boundaryMesh_(bm)
99{
100 if (constraintType(patchType))
101 {
102 addGroup(patchType);
103 }
104}
105
106
108(
109 const word& name,
110 const label size,
111 const label start,
112 const label index,
113 const polyBoundaryMesh& bm,
114 const word& physicalType,
115 const wordList& inGroups
116)
117:
118 patchIdentifier(name, index, physicalType, inGroups),
120 (
121 SubList<face>(bm.mesh().faces(), size, start),
123 ),
124 start_(start),
125 boundaryMesh_(bm)
126{}
127
128
130(
131 const word& name,
132 const dictionary& dict,
133 const label index,
134 const polyBoundaryMesh& bm,
135 const word& patchType
136)
137:
138 patchIdentifier(name, dict, index),
140 (
142 (
143 bm.mesh().faces(),
144 dict.get<label>("nFaces"),
145 dict.get<label>("startFace")
146 ),
147 bm.mesh().points()
148 ),
149 start_(dict.get<label>("startFace")),
150 boundaryMesh_(bm)
151{
152 if (constraintType(patchType))
153 {
154 addGroup(patchType);
155 }
156}
157
158
160(
161 const polyPatch& pp,
162 const polyBoundaryMesh& bm
163)
164:
165 patchIdentifier(pp),
167 (
168 SubList<face>
169 (
170 bm.mesh().faces(),
171 pp.size(),
172 pp.start()
173 ),
175 ),
176 start_(pp.start()),
177 boundaryMesh_(bm)
178{}
179
180
182(
183 const polyPatch& pp,
184 const polyBoundaryMesh& bm,
185 const label index,
186 const label newSize,
187 const label newStart
188)
189:
190 patchIdentifier(pp, index),
192 (
194 (
195 bm.mesh().faces(),
196 newSize,
197 newStart
198 ),
200 ),
201 start_(newStart),
202 boundaryMesh_(bm)
203{}
204
205
207(
208 const polyPatch& pp,
209 const polyBoundaryMesh& bm,
210 const label index,
211 const labelUList& mapAddressing,
212 const label newStart
213)
214:
215 patchIdentifier(pp, index),
217 (
219 (
220 bm.mesh().faces(),
221 mapAddressing.size(),
222 newStart
223 ),
225 ),
226 start_(newStart),
227 boundaryMesh_(bm)
228{}
229
230
232:
235 start_(p.start_),
236 boundaryMesh_(p.boundaryMesh_)
237{}
238
239
241(
242 const polyPatch& p,
243 const labelList& faceCells
244)
245:
246 polyPatch(p)
248 faceCellsPtr_.reset(new labelList::subList(faceCells, faceCells.size()));
249}
250
251
252// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
253
257}
258
259
260// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
261
262bool Foam::polyPatch::constraintType(const word& patchType)
263{
264 return
265 (
266 !patchType.empty()
269 );
270}
271
272
274{
275 const auto& cnstrTable = *dictionaryConstructorTablePtr_;
276
277 wordList cTypes(cnstrTable.size());
278
279 label i = 0;
280
281 forAllConstIters(cnstrTable, iter)
282 {
283 if (constraintType(iter.key()))
284 {
285 cTypes[i++] = iter.key();
286 }
287 }
289 cTypes.resize(i);
290
291 return cTypes;
292}
293
294
296{
297 // Same as start_ - polyMesh::nInternalFaces()
298 return start_ - boundaryMesh_.start();
299}
300
303{
304 return boundaryMesh_;
305}
306
311}
312
313
315{
316 return patchSlice(boundaryMesh().mesh().faceOwner());
317}
318
319
320// Potentially useful to simplify logic elsewhere?
321// const Foam::labelList::subList Foam::polyPatch::faceNeighbour() const
322// {
323// return labelList::subList();
324// }
325
330}
331
334{
336}
337
338
340{
341 auto tcc = tmp<vectorField>::New(size());
342 auto& cc = tcc.ref();
343
344 // get reference to global cell centres
345 const vectorField& gcc = boundaryMesh_.mesh().cellCentres();
346
347 const labelUList& faceCells = this->faceCells();
348
349 forAll(faceCells, facei)
350 {
351 cc[facei] = gcc[faceCells[facei]];
352 }
353
354 return tcc;
355}
356
357
359(
360 const pointField& points
361) const
362{
363 auto tfraction = tmp<scalarField>::New(size());
364 auto& fraction = tfraction.ref();
365
366 const vectorField::subField faceAreas = this->faceAreas();
367
368 forAll(*this, facei)
369 {
370 const face& f = this->operator[](facei);
371 fraction[facei] = faceAreas[facei].mag()/(f.mag(points) + ROOTVSMALL);
372 }
373
374 return tfraction;
375}
376
377
379{
380 if (areaFractionPtr_)
382 return tmp<scalarField>(*areaFractionPtr_);
383 }
384 return nullptr;
385}
386
388void Foam::polyPatch::areaFraction(const scalar fraction)
389{
390 areaFractionPtr_ = std::make_unique<scalarField>(size(), fraction);
391}
392
393
395{
396 if (fraction)
397 {
398 // Steal or clone
399 areaFractionPtr_.reset(fraction.ptr());
400 }
401 else
402 {
403 areaFractionPtr_.reset(nullptr);
404 }
405}
406
407
409{
410 if (!faceCellsPtr_)
411 {
412 faceCellsPtr_.reset
413 (
415 (
416 patchSlice(boundaryMesh().mesh().faceOwner())
417 )
418 );
419 }
420
421 return *faceCellsPtr_;
422}
423
424
426{
427 if (!mePtr_)
428 {
429 mePtr_.reset
430 (
431 new labelList
432 (
434 (
435 boundaryMesh().mesh().edges(),
436 boundaryMesh().mesh().pointEdges()
437 )
438 )
439 );
440 }
441
442 return *mePtr_;
443}
444
445
447{
450 faceCellsPtr_.reset(nullptr);
451 mePtr_.reset(nullptr);
452 areaFractionPtr_.reset(nullptr);
453}
454
455
457{
458 os.writeEntry("type", type());
460 os.writeEntry("nFaces", size());
461 os.writeEntry("startFace", start());
463
464
466{}
467
468
470(
471 PstreamBuffers&,
472 const primitivePatch&,
474 labelList& rotation
475) const
476{
477 // Nothing changed.
478 return false;
479}
480
481
482// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
483
485{
486 clearAddressing();
487
490 start_ = p.start_;
491}
492
493
494// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
495
496Foam::Ostream& Foam::operator<<(Ostream& os, const polyPatch& p)
497{
498 p.write(os);
499 os.check(FUNCTION_NAME);
500 return os;
501}
502
503
504// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
SubField< vector > subField
Definition Field.H:183
SubList< label > subList
Definition List.H:129
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
void operator=(const PrimitivePatch< SubList< face >, const pointField & > &rhs)
const Field< point_type > & points() const noexcept
virtual void movePoints(const Field< point_type > &)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
UList< T > & reset(std::nullptr_t) noexcept
Reset to zero-sized and nullptr.
Definition SubListI.H:109
bool get(const label i) const
Definition UList.H:868
void size(const label n)
Definition UList.H:118
face & operator[](const label i)
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
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Identifies a patch by name and index, with optional physical type and group information.
const word & physicalType() const noexcept
The (optional) physical type of the patch.
patchIdentifier & operator=(const patchIdentifier &)=default
Copy assignment.
const wordList & inGroups() const noexcept
The (optional) groups that the patch belongs to.
void write(Ostream &os) const
Write (physicalType, inGroups) dictionary entries (without surrounding braces).
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
patchIdentifier(const patchIdentifier &)=default
Copy construct.
void addGroup(const word &name)
Add (unique) group for the patch.
Abstract base class for point-mesh patch fields.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition polyPatch.C:266
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition polyPatch.C:59
const vectorField::subField faceAreas() const
Return face normals.
Definition polyPatch.C:326
friend class polyBoundaryMesh
Definition polyPatch.H:112
virtual void clearGeom()
Clear geometry.
Definition polyPatch.C:66
void operator=(const polyPatch &p)
Copy assignment.
Definition polyPatch.C:477
label offset() const noexcept
The offset where this patch starts in the boundary face list.
Definition polyPatch.C:288
const faceList::subList faces() const
Return mesh faces for the patch.
Definition polyPatch.C:301
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
Definition polyPatch.C:458
virtual void write(Ostream &os) const
Write the polyPatch data as a dictionary.
Definition polyPatch.C:449
const vectorField::subField faceCentres() const
Return face centres.
Definition polyPatch.C:320
virtual void movePoints(PstreamBuffers &, const pointField &p)
Correct patches after moving points.
Definition polyPatch.C:53
virtual ~polyPatch()
Destructor.
Definition polyPatch.C:247
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition polyPatch.C:295
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition polyPatch.C:463
virtual void clearAddressing()
Clear addressing.
Definition polyPatch.C:439
polyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType)
Construct from components.
Definition polyPatch.C:75
static int disallowGenericPolyPatch
Debug switch to disallow the use of genericPolyPatch.
Definition polyPatch.H:164
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
const List< T >::subList patchSlice(const UList< T > &values) const
This patch slice from the complete list, which has size mesh::nFaces(), using the number of patch fac...
Definition polyPatch.H:502
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition polyPatch.C:332
const labelList::subList faceOwner() const
Return face owner for the patch.
Definition polyPatch.C:307
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition polyPatch.C:255
const labelUList & faceCells() const
Return face-cell addressing.
Definition polyPatch.C:401
const labelList & meshEdges() const
Return global edge index for local edges.
Definition polyPatch.C:418
tmp< scalarField > areaFraction() const
Return the cached area fraction. Usually only set for the non-overlap patches on ACMI.
Definition polyPatch.C:371
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition tmpI.H:256
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
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
auto & name
const pointField & points
#define FUNCTION_NAME
int debugSwitch(const char *name, const int deflt=0)
Lookup debug switch or add default value.
Definition debug.C:222
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Field< vector > vectorField
Specialisation of Field<T> for vector.
const direction noexcept
Definition scalarImpl.H:265
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
labelList f(nPoints)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235