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-2022 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 #include "demandDrivenData.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(polyPatch, 0);
45 
47  (
48  debug::debugSwitch("disallowGenericPolyPatch", 0)
49  );
50 
51  defineRunTimeSelectionTable(polyPatch, word);
52  defineRunTimeSelectionTable(polyPatch, dictionary);
53 
56 }
57 
58 
59 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
60 
62 {
64 }
65 
66 
68 {
70  clearAddressing();
71 }
72 
73 
75 {
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const word& name,
85  const label size,
86  const label start,
87  const label index,
88  const polyBoundaryMesh& bm,
89  const word& patchType
90 )
91 :
92  patchIdentifier(name, index),
94  (
95  faceSubList(bm.mesh().faces(), size, start),
96  bm.mesh().points()
97  ),
98  start_(start),
99  boundaryMesh_(bm),
100  faceCellsPtr_(nullptr),
101  mePtr_(nullptr)
102 {
103  if (constraintType(patchType))
104  {
105  inGroups().appendUniq(patchType);
106  }
107 }
108 
109 
111 (
112  const word& name,
113  const label size,
114  const label start,
115  const label index,
116  const polyBoundaryMesh& bm,
117  const word& physicalType,
118  const wordList& inGroups
119 )
120 :
121  patchIdentifier(name, index, physicalType, inGroups),
123  (
124  faceSubList(bm.mesh().faces(), size, start),
125  bm.mesh().points()
126  ),
127  start_(start),
128  boundaryMesh_(bm),
129  faceCellsPtr_(nullptr),
130  mePtr_(nullptr)
131 {}
132 
133 
135 (
136  const word& name,
137  const dictionary& dict,
138  const label index,
139  const polyBoundaryMesh& bm,
140  const word& patchType
141 )
142 :
143  patchIdentifier(name, dict, index),
145  (
147  (
148  bm.mesh().faces(),
149  dict.get<label>("nFaces"),
150  dict.get<label>("startFace")
151  ),
152  bm.mesh().points()
153  ),
154  start_(dict.get<label>("startFace")),
155  boundaryMesh_(bm),
156  faceCellsPtr_(nullptr),
157  mePtr_(nullptr)
158 {
159  if (constraintType(patchType))
160  {
161  inGroups().appendUniq(patchType);
162  }
163 }
164 
165 
167 (
168  const polyPatch& pp,
169  const polyBoundaryMesh& bm
170 )
171 :
172  patchIdentifier(pp),
174  (
176  (
177  bm.mesh().faces(),
178  pp.size(),
179  pp.start()
180  ),
181  bm.mesh().points()
182  ),
183  start_(pp.start()),
184  boundaryMesh_(bm),
185  faceCellsPtr_(nullptr),
186  mePtr_(nullptr)
187 {}
188 
189 
191 (
192  const polyPatch& pp,
193  const polyBoundaryMesh& bm,
194  const label index,
195  const label newSize,
196  const label newStart
197 )
198 :
199  patchIdentifier(pp, index),
201  (
203  (
204  bm.mesh().faces(),
205  newSize,
206  newStart
207  ),
208  bm.mesh().points()
209  ),
210  start_(newStart),
211  boundaryMesh_(bm),
212  faceCellsPtr_(nullptr),
213  mePtr_(nullptr)
214 {}
215 
216 
218 (
219  const polyPatch& pp,
220  const polyBoundaryMesh& bm,
221  const label index,
222  const labelUList& mapAddressing,
223  const label newStart
224 )
225 :
226  patchIdentifier(pp, index),
228  (
230  (
231  bm.mesh().faces(),
232  mapAddressing.size(),
233  newStart
234  ),
235  bm.mesh().points()
236  ),
237  start_(newStart),
238  boundaryMesh_(bm),
239  faceCellsPtr_(nullptr),
240  mePtr_(nullptr)
241 {}
242 
243 
245 :
247  primitivePatch(p),
248  start_(p.start_),
249  boundaryMesh_(p.boundaryMesh_),
250  faceCellsPtr_(nullptr),
251  mePtr_(nullptr)
252 {}
253 
254 
256 (
257  const polyPatch& p,
258  const labelList& faceCells
259 )
260 :
261  polyPatch(p)
262 {
263  faceCellsPtr_ = new labelList::subList(faceCells, faceCells.size());
264 }
265 
266 
267 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
268 
270 {
271  clearAddressing();
272 }
273 
274 
275 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
276 
277 bool Foam::polyPatch::constraintType(const word& patchType)
278 {
279  return
280  (
281  !patchType.empty()
284  );
285 }
286 
287 
289 {
290  const auto& cnstrTable = *dictionaryConstructorTablePtr_;
291 
292  wordList cTypes(cnstrTable.size());
293 
294  label i = 0;
295 
296  forAllConstIters(cnstrTable, iter)
297  {
298  if (constraintType(iter.key()))
299  {
300  cTypes[i++] = iter.key();
301  }
302  }
304  cTypes.resize(i);
305 
306  return cTypes;
307 }
308 
310 Foam::label Foam::polyPatch::offset() const
311 {
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_ = new labelList::subList
377  (
378  patchSlice(boundaryMesh().mesh().faceOwner())
379  );
380  }
381 
382  return *faceCellsPtr_;
383 }
384 
385 
387 {
388  if (!mePtr_)
389  {
390  mePtr_ =
391  new labelList
392  (
394  (
395  boundaryMesh().mesh().edges(),
396  boundaryMesh().mesh().pointEdges()
397  )
398  );
399  }
400 
401  return *mePtr_;
402 }
403 
404 
406 {
409  deleteDemandDrivenData(faceCellsPtr_);
410  deleteDemandDrivenData(mePtr_);
411 }
412 
413 
414 void Foam::polyPatch::write(Ostream& os) const
415 {
416  os.writeEntry("type", type());
418  os.writeEntry("nFaces", size());
419  os.writeEntry("startFace", start());
420 }
422 
424 {}
425 
426 
428 (
430  const primitivePatch&,
432  labelList& rotation
433 ) const
434 {
435  // Nothing changed.
436  return false;
437 }
438 
439 
440 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
441 
443 {
444  clearAddressing();
445 
448  start_ = p.start_;
449 }
450 
451 
452 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
453 
454 Foam::Ostream& Foam::operator<<(Ostream& os, const polyPatch& p)
455 {
456  p.write(os);
458  return os;
459 }
460 
461 
462 // ************************************************************************* //
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:120
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.
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:309
const wordList & inGroups() const noexcept
The (optional) groups that the patch belongs to.
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
label appendUniq(const T &val)
Append an element if not already in the list.
Definition: List.H:526
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
SubList< label > subList
Declare type of subList.
Definition: List.H:122
virtual void clearAddressing()
Clear addressing.
Definition: polyPatch.C:398
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:62
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:67
void operator=(const PrimitivePatch< FaceList, PointField > &rhs)
Copy assign faces. Leave points alone (could be a reference).
SubList< face > faceSubList
A SubList of faces.
Definition: faceListFwd.H:42
virtual void movePoints(PstreamBuffers &, const pointField &p)
Correct patches after moving points.
Definition: polyPatch.C:54
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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:752
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:50
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 INVALID.
Definition: exprTraits.C:52
virtual void write(Ostream &os) const
Write the polyPatch data as a dictionary.
Definition: polyPatch.C:407
const pointField & points
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:262
label offset() const
The offset where this patch starts in the boundary face list.
Definition: polyPatch.C:303
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:327
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:60
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:76
void operator=(const polyPatch &)
Assignment.
Definition: polyPatch.C:435
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:55
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:379
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...
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: polyPatch.C:281
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:76
List< word > wordList
A List of words.
Definition: fileName.H:58
Template functions to aid in the implementation of demand driven data.
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:270
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: polyPatch.C:421
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.
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:416
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
virtual bool write(const bool valid=true) const
Write using setting from DB.
void deleteDemandDrivenData(DataPtr &dataPtr)
patchIdentifier & operator=(const patchIdentifier &)=default
Copy assignment.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:315