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-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(polyPatch, 0);
44 
46  (
47  debug::debugSwitch("disallowGenericPolyPatch", 0)
48  );
49 
50  defineRunTimeSelectionTable(polyPatch, word);
51  defineRunTimeSelectionTable(polyPatch, dictionary);
52 
55 }
56 
57 
58 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
59 
61 {
63 }
64 
65 
67 {
69  clearAddressing();
70 }
71 
72 
74 {
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  faceSubList(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  faceSubList(bm.mesh().faces(), size, start),
122  bm.mesh().points()
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  (
169  (
170  bm.mesh().faces(),
171  pp.size(),
172  pp.start()
173  ),
174  bm.mesh().points()
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  ),
199  bm.mesh().points()
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  ),
224  bm.mesh().points()
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)
247 {
248  faceCellsPtr_.reset(new labelList::subList(faceCells, faceCells.size()));
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
253 
255 {
256  clearAddressing();
257 }
258 
259 
260 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
261 
262 bool 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 
309 {
310  return patchSlice(boundaryMesh().mesh().faceCentres());
311 }
312 
315 {
316  return patchSlice(boundaryMesh().mesh().faceAreas());
317 }
318 
319 
321 {
322  auto tcc = tmp<vectorField>::New(size());
323  auto& cc = tcc.ref();
324 
325  // get reference to global cell centres
326  const vectorField& gcc = boundaryMesh_.mesh().cellCentres();
327 
328  const labelUList& faceCells = this->faceCells();
329 
330  forAll(faceCells, facei)
331  {
332  cc[facei] = gcc[faceCells[facei]];
333  }
334 
335  return tcc;
336 }
337 
338 
340 (
341  const pointField& points
342 ) const
343 {
344  auto tfraction = tmp<scalarField>::New(size());
345  auto& fraction = tfraction.ref();
346 
347  const vectorField::subField faceAreas = this->faceAreas();
348 
349  forAll(*this, facei)
350  {
351  const face& f = this->operator[](facei);
352  fraction[facei] = faceAreas[facei].mag()/(f.mag(points) + ROOTVSMALL);
353  }
354 
355  return tfraction;
356 }
357 
358 
360 {
361  if (areaFractionPtr_)
362  {
363  return tmp<scalarField>(*areaFractionPtr_);
364  }
365  return nullptr;
366 }
367 
369 void Foam::polyPatch::areaFraction(const scalar fraction)
370 {
371  areaFractionPtr_ = std::make_unique<scalarField>(size(), fraction);
372 }
373 
374 
376 {
377  if (fraction)
378  {
379  // Steal or clone
380  areaFractionPtr_.reset(fraction.ptr());
381  }
382  else
383  {
384  areaFractionPtr_.reset(nullptr);
385  }
386 }
387 
388 
390 {
391  if (!faceCellsPtr_)
392  {
393  faceCellsPtr_.reset
394  (
396  (
397  patchSlice(boundaryMesh().mesh().faceOwner())
398  )
399  );
400  }
401 
402  return *faceCellsPtr_;
403 }
404 
405 
407 {
408  if (!mePtr_)
409  {
410  mePtr_.reset
411  (
412  new labelList
413  (
415  (
416  boundaryMesh().mesh().edges(),
417  boundaryMesh().mesh().pointEdges()
418  )
419  )
420  );
421  }
422 
423  return *mePtr_;
424 }
425 
426 
428 {
431  faceCellsPtr_.reset(nullptr);
432  mePtr_.reset(nullptr);
433  areaFractionPtr_.reset(nullptr);
434 }
435 
436 
437 void Foam::polyPatch::write(Ostream& os) const
438 {
439  os.writeEntry("type", type());
441  os.writeEntry("nFaces", size());
442  os.writeEntry("startFace", start());
443 }
445 
447 {}
448 
449 
451 (
453  const primitivePatch&,
455  labelList& rotation
456 ) const
457 {
458  // Nothing changed.
459  return false;
460 }
461 
462 
463 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
464 
466 {
467  clearAddressing();
468 
471  start_ = p.start_;
472 }
473 
474 
475 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
476 
477 Foam::Ostream& Foam::operator<<(Ostream& os, const polyPatch& p)
478 {
479  p.write(os);
481  return os;
482 }
483 
484 
485 // ************************************************************************* //
void addGroup(const word &name)
Add (unique) group for the patch.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
void write(Ostream &os) const
Write (physicalType, inGroups) dictionary entries (without surrounding braces)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static int disallowGenericPolyPatch
Debug switch to disallow the use of genericPolyPatch.
Definition: polyPatch.H:164
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
Identifies a patch by name and index, with optional physical type and group information.
int debugSwitch(const char *name, const int deflt=0)
Lookup debug switch or add default value.
Definition: debug.C:222
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
SubList< label > subList
Declare type of subList.
Definition: List.H:144
virtual void clearAddressing()
Clear addressing.
Definition: polyPatch.C:420
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
SubField is a Field obtained as a section of another Field, without its own allocation. SubField is derived from a SubList rather than a List.
Definition: Field.H:63
Macros for easy insertion into run-time selection tables.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual void clearGeom()
Clear geometry.
Definition: polyPatch.C:66
void operator=(const PrimitivePatch< FaceList, PointField > &rhs)
Copy assign faces. Leave points alone (could be a reference).
virtual void movePoints(PstreamBuffers &, const pointField &p)
Correct patches after moving points.
Definition: polyPatch.C:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< scalarField > areaFraction() const
Return the cached area fraction. Usually only set for the non-overlap patches on ACMI.
Definition: polyPatch.C:352
Abstract base class for point-mesh patch fields.
virtual void movePoints(const Field< point_type > &)
Correct patch after moving points.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:50
void operator=(const polyPatch &p)
Copy assignment.
Definition: polyPatch.C:458
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual void write(Ostream &os) const
Write the polyPatch data as a dictionary.
Definition: polyPatch.C:430
const pointField & points
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: polyPatch.C:295
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:382
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual ~polyPatch()
Destructor.
Definition: polyPatch.C:247
label offset() const noexcept
The offset where this patch starts in the boundary face list.
Definition: polyPatch.C:288
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:313
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:59
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
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:307
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
defineTypeNameAndDebug(combustionModel, 0)
const labelList & meshEdges() const
Return global edge index for local edges.
Definition: polyPatch.C:399
labelList f(nPoints)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: polyPatch.C:266
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
List< word > wordList
List of word.
Definition: fileName.H:59
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:255
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: polyPatch.C:444
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using cell addressing.
SubList< face > faceSubList
SubList of faces.
Definition: faceListFwd.H:42
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
Definition: polyPatch.C:439
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:256
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
patchIdentifier & operator=(const patchIdentifier &)=default
Copy assignment.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:301