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-2023 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  faceCellsPtr_(nullptr),
100  mePtr_(nullptr)
101 {
102  if (constraintType(patchType))
103  {
104  addGroup(patchType);
105  }
106 }
107 
108 
110 (
111  const word& name,
112  const label size,
113  const label start,
114  const label index,
115  const polyBoundaryMesh& bm,
116  const word& physicalType,
117  const wordList& inGroups
118 )
119 :
120  patchIdentifier(name, index, physicalType, inGroups),
122  (
123  faceSubList(bm.mesh().faces(), size, start),
124  bm.mesh().points()
125  ),
126  start_(start),
127  boundaryMesh_(bm),
128  faceCellsPtr_(nullptr),
129  mePtr_(nullptr)
130 {}
131 
132 
134 (
135  const word& name,
136  const dictionary& dict,
137  const label index,
138  const polyBoundaryMesh& bm,
139  const word& patchType
140 )
141 :
142  patchIdentifier(name, dict, index),
144  (
146  (
147  bm.mesh().faces(),
148  dict.get<label>("nFaces"),
149  dict.get<label>("startFace")
150  ),
151  bm.mesh().points()
152  ),
153  start_(dict.get<label>("startFace")),
154  boundaryMesh_(bm),
155  faceCellsPtr_(nullptr),
156  mePtr_(nullptr)
157 {
158  if (constraintType(patchType))
159  {
160  addGroup(patchType);
161  }
162 }
163 
164 
166 (
167  const polyPatch& pp,
168  const polyBoundaryMesh& bm
169 )
170 :
171  patchIdentifier(pp),
173  (
175  (
176  bm.mesh().faces(),
177  pp.size(),
178  pp.start()
179  ),
180  bm.mesh().points()
181  ),
182  start_(pp.start()),
183  boundaryMesh_(bm),
184  faceCellsPtr_(nullptr),
185  mePtr_(nullptr)
186 {}
187 
188 
190 (
191  const polyPatch& pp,
192  const polyBoundaryMesh& bm,
193  const label index,
194  const label newSize,
195  const label newStart
196 )
197 :
198  patchIdentifier(pp, index),
200  (
202  (
203  bm.mesh().faces(),
204  newSize,
205  newStart
206  ),
207  bm.mesh().points()
208  ),
209  start_(newStart),
210  boundaryMesh_(bm),
211  faceCellsPtr_(nullptr),
212  mePtr_(nullptr)
213 {}
214 
215 
217 (
218  const polyPatch& pp,
219  const polyBoundaryMesh& bm,
220  const label index,
221  const labelUList& mapAddressing,
222  const label newStart
223 )
224 :
225  patchIdentifier(pp, index),
227  (
229  (
230  bm.mesh().faces(),
231  mapAddressing.size(),
232  newStart
233  ),
234  bm.mesh().points()
235  ),
236  start_(newStart),
237  boundaryMesh_(bm),
238  faceCellsPtr_(nullptr),
239  mePtr_(nullptr)
240 {}
241 
242 
244 :
246  primitivePatch(p),
247  start_(p.start_),
248  boundaryMesh_(p.boundaryMesh_),
249  faceCellsPtr_(nullptr),
250  mePtr_(nullptr)
251 {}
252 
253 
255 (
256  const polyPatch& p,
257  const labelList& faceCells
258 )
259 :
260  polyPatch(p)
261 {
262  faceCellsPtr_.reset(new labelList::subList(faceCells, faceCells.size()));
263 }
264 
265 
266 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
267 
269 {
270  clearAddressing();
271 }
272 
273 
274 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
275 
276 bool Foam::polyPatch::constraintType(const word& patchType)
277 {
278  return
279  (
280  !patchType.empty()
283  );
284 }
285 
286 
288 {
289  const auto& cnstrTable = *dictionaryConstructorTablePtr_;
290 
291  wordList cTypes(cnstrTable.size());
292 
293  label i = 0;
294 
295  forAllConstIters(cnstrTable, iter)
296  {
297  if (constraintType(iter.key()))
298  {
299  cTypes[i++] = iter.key();
300  }
301  }
303  cTypes.resize(i);
304 
305  return cTypes;
306 }
307 
308 
310 {
311  // Same as start_ - polyMesh::nInternalFaces()
312  return start_ - boundaryMesh_.start();
313 }
314 
317 {
318  return boundaryMesh_;
319 }
320 
323 {
324  return patchSlice(boundaryMesh().mesh().faceCentres());
325 }
326 
329 {
330  return patchSlice(boundaryMesh().mesh().faceAreas());
331 }
332 
333 
335 {
336  tmp<vectorField> tcc(new vectorField(size()));
337  vectorField& cc = tcc.ref();
338 
339  // get reference to global cell centres
340  const vectorField& gcc = boundaryMesh_.mesh().cellCentres();
341 
342  const labelUList& faceCells = this->faceCells();
343 
344  forAll(faceCells, facei)
345  {
346  cc[facei] = gcc[faceCells[facei]];
347  }
348 
349  return tcc;
350 }
351 
352 
354 {
355  tmp<scalarField> tfraction(new scalarField(size()));
356  scalarField& fraction = tfraction.ref();
357 
358  const vectorField::subField faceAreas = this->faceAreas();
359  const pointField& points = this->points();
360 
361  forAll(*this, facei)
362  {
363  const face& curFace = this->operator[](facei);
364  fraction[facei] =
365  mag(faceAreas[facei])/(curFace.mag(points) + ROOTVSMALL);
366  }
367 
368  return tfraction;
369 }
370 
371 
373 {
374  if (!faceCellsPtr_)
375  {
376  faceCellsPtr_.reset
377  (
379  (
380  patchSlice(boundaryMesh().mesh().faceOwner())
381  )
382  );
383  }
384 
385  return *faceCellsPtr_;
386 }
387 
388 
390 {
391  if (!mePtr_)
392  {
393  mePtr_.reset
394  (
395  new labelList
396  (
398  (
399  boundaryMesh().mesh().edges(),
400  boundaryMesh().mesh().pointEdges()
401  )
402  )
403  );
404  }
405 
406  return *mePtr_;
407 }
408 
409 
411 {
414  faceCellsPtr_.reset(nullptr);
415  mePtr_.reset(nullptr);
416 }
417 
418 
419 void Foam::polyPatch::write(Ostream& os) const
420 {
421  os.writeEntry("type", type());
423  os.writeEntry("nFaces", size());
424  os.writeEntry("startFace", start());
425 }
427 
429 {}
430 
431 
433 (
435  const primitivePatch&,
437  labelList& rotation
438 ) const
439 {
440  // Nothing changed.
441  return false;
442 }
443 
444 
445 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
446 
448 {
449  clearAddressing();
450 
453  start_ = p.start_;
454 }
455 
456 
457 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
458 
459 Foam::Ostream& Foam::operator<<(Ostream& os, const polyPatch& p)
460 {
461  p.write(os);
463  return os;
464 }
465 
466 
467 // ************************************************************************* //
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)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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:159
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:403
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 area fraction as the ratio of the stored face area and the area given by the face points...
Definition: polyPatch.C:346
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:440
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:412
const pointField & points
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: polyPatch.C:309
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:365
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual ~polyPatch()
Destructor.
Definition: polyPatch.C:261
label offset() const noexcept
The offset where this patch starts in the boundary face list.
Definition: polyPatch.C:302
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:327
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:321
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:382
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:280
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:269
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: polyPatch.C:426
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
Definition: polyPatch.C:421
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:315