polyPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "polyPatch.H"
28 #include "polyBoundaryMesh.H"
29 #include "polyMesh.H"
30 #include "primitiveMesh.H"
31 #include "SubField.H"
32 #include "entry.H"
33 #include "dictionary.H"
34 #include "pointPatchField.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(polyPatch, 0);
41 
43  (
44  debug::debugSwitch("disallowGenericPolyPatch", 0)
45  );
46 
47  defineRunTimeSelectionTable(polyPatch, word);
48  defineRunTimeSelectionTable(polyPatch, dictionary);
49 
50  addToRunTimeSelectionTable(polyPatch, polyPatch, word);
51  addToRunTimeSelectionTable(polyPatch, polyPatch, dictionary);
52 }
53 
54 
55 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
56 
58 {
60 }
61 
62 
64 {
66 }
67 
68 
70 {
72  clearAddressing();
73 }
74 
75 
77 {
79 }
80 
81 
82 void Foam::polyPatch::rename(const wordList& newNames)
83 {
84  name_ = newNames[index_];
85 }
86 
87 
88 void Foam::polyPatch::reorder(const labelUList& newToOldIndex)
89 {
90  index_ = findIndex(newToOldIndex, index_);
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const word& name,
99  const label size,
100  const label start,
101  const label index,
102  const polyBoundaryMesh& bm,
103  const word& patchType
104 )
105 :
106  patchIdentifier(name, index),
108  (
109  faceSubList(bm.mesh().faces(), size, start),
110  bm.mesh().points()
111  ),
112  start_(start),
113  boundaryMesh_(bm),
114  faceCellsPtr_(nullptr),
115  mePtr_(nullptr)
116 {
117  if
118  (
119  patchType != word::null
120  && constraintType(patchType)
121  && findIndex(inGroups(), patchType) == -1
122  )
123  {
124  inGroups().append(patchType);
125  }
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.lookup<label>("nFaces"),
145  dict.lookup<label>("startFace")
146  ),
147  bm.mesh().points()
148  ),
149  start_(dict.lookup<label>("startFace")),
150  boundaryMesh_(bm),
151  faceCellsPtr_(nullptr),
152  mePtr_(nullptr)
153 {
154  if
155  (
156  patchType != word::null
157  && constraintType(patchType)
158  && findIndex(inGroups(), patchType) == -1
159  )
160  {
161  inGroups().append(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 :
246  patchIdentifier(p),
247  primitivePatch(p),
248  start_(p.start_),
249  boundaryMesh_(p.boundaryMesh_),
250  faceCellsPtr_(nullptr),
251  mePtr_(nullptr)
252 {}
253 
254 
255 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
256 
258 {
259  clearAddressing();
260 }
261 
262 
263 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
264 
266 {
268 }
269 
270 
272 {
273  wordList cTypes(dictionaryConstructorTablePtr_->size());
274 
275  label i = 0;
276 
277  for
278  (
279  dictionaryConstructorTable::iterator cstrIter =
280  dictionaryConstructorTablePtr_->begin();
281  cstrIter != dictionaryConstructorTablePtr_->end();
282  ++cstrIter
283  )
284  {
285  if (constraintType(cstrIter.key()))
286  {
287  cTypes[i++] = cstrIter.key();
288  }
289  }
290 
291  cTypes.setSize(i);
292 
293  return cTypes;
294 }
295 
296 
298 {
299  return boundaryMesh_;
300 }
301 
302 
304 {
305  return patchSlice(boundaryMesh().mesh().faceCentres());
306 }
307 
308 
310 {
311  return patchSlice(boundaryMesh().mesh().faceAreas());
312 }
313 
314 
316 {
318 }
319 
320 
322 {
323  tmp<vectorField> tcc(new vectorField(size()));
324  vectorField& cc = tcc.ref();
325 
326  // get reference to global cell centres
327  const vectorField& gcc = boundaryMesh_.mesh().cellCentres();
328 
329  const labelUList& faceCells = this->faceCells();
330 
331  forAll(faceCells, facei)
332  {
333  cc[facei] = gcc[faceCells[facei]];
334  }
335 
336  return tcc;
337 }
338 
339 
341 {
342  if (!faceCellsPtr_)
343  {
344  faceCellsPtr_ = new labelList::subList
345  (
346  patchSlice(boundaryMesh().mesh().faceOwner())
347  );
348  }
349 
350  return *faceCellsPtr_;
351 }
352 
353 
355 {
356  if (!mePtr_)
357  {
358  mePtr_ =
359  new labelList
360  (
362  (
363  boundaryMesh().mesh().edges(),
365  )
366  );
367  }
368 
369  return *mePtr_;
370 }
371 
372 
374 {
377  deleteDemandDrivenData(faceCellsPtr_);
378  deleteDemandDrivenData(mePtr_);
379 }
380 
381 
383 {
384  writeEntry(os, "type", type());
386  writeEntry(os, "nFaces", size());
387  writeEntry(os, "startFace", start());
388 }
389 
390 
392 {}
393 
394 
396 (
398  const primitivePatch&,
399  labelList& faceMap,
400  labelList& rotation
401 ) const
402 {
403  // Nothing changed.
404  return false;
405 }
406 
407 
408 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
409 
411 {
412  clearAddressing();
413 
414  patchIdentifier::operator=(p);
416  start_ = p.start_;
417 }
418 
419 
420 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
421 
423 {
424  p.write(os);
425  os.check("Ostream& operator<<(Ostream& os, const polyPatch& p");
426  return os;
427 }
428 
429 
430 // ************************************************************************* //
const labelListList & pointEdges() const
Return point-edge addressing.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static int disallowGenericPolyPatch
Debug switch to disallow the use of genericPolyPatch.
Definition: polyPatch.H:137
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Identifies patch by name, patch index and physical type.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:297
SubList< label > subList
Declare type of subList.
Definition: List.H:199
virtual void clearAddressing()
Clear addressing.
Definition: polyPatch.C:373
Pre-declare related SubField type.
Definition: Field.H:60
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
virtual void clearGeom()
Clear geometry.
Definition: polyPatch.C:76
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: polyPatch.C:382
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:265
Abstract base class for point-mesh patch fields.
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:53
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
A class for handling words, derived from string.
Definition: word.H:59
virtual ~polyPatch()
Destructor.
Definition: polyPatch.C:257
virtual void rename(const wordList &newNames)
Reset the patch name.
Definition: polyPatch.C:82
int debugSwitch(const char *name, const int defaultValue=0)
Lookup debug switch or add default value.
Definition: debug.C:211
static const word null
An empty word.
Definition: word.H:77
const polyMesh & mesh() const
Return the mesh reference.
virtual void reorder(const labelUList &newToOldIndex)
Reset the patch index.
Definition: polyPatch.C:88
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:321
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
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:97
const vectorField & cellCentres() const
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
void operator=(const polyPatch &)
Assignment.
Definition: polyPatch.C:410
Foam::polyBoundaryMesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const vectorField::subField faceAreas() const
Return face areas.
Definition: polyPatch.C:309
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
const labelList & meshEdges() const
Return global edge index for local edges.
Definition: polyPatch.C:354
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
virtual void topoChange(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:69
void write(Ostream &) const
Write patchIdentifier as a dictionary.
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: polyPatch.C:271
virtual void movePoints(const pointField &p)
Correct patches after moving points.
Definition: polyPatch.C:57
const scalarField::subField magFaceAreas() const
Return face area magnitudes.
Definition: polyPatch.C:315
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: polyPatch.C:396
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:306
Ostream & operator<<(Ostream &, const ensightPart &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SubList< face > faceSubList
Definition: faceListFwd.H:44
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialise ordering for primitivePatch. Does not.
Definition: polyPatch.C:391
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
void deleteDemandDrivenData(DataPtr &dataPtr)
void operator=(const PrimitivePatch< FaceList, PointField > &)
Assignment operator.
const List< T >::subList patchSlice(const UList< T > &l) const
Slice list to patch.
Definition: polyPatch.H:339
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:303