faPatch.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2019-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 "faPatch.H"
31 #include "faBoundaryMesh.H"
32 #include "faMesh.H"
33 #include "areaFields.H"
34 #include "edgeFields.H"
35 #include "edgeHashes.H"
36 #include "polyMesh.H"
37 #include "polyPatch.H"
38 //#include "pointPatchField.H"
39 #include "demandDrivenData.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(faPatch, 0);
48 }
49 
50 
51 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
52 
53 bool Foam::faPatch::constraintType(const word& patchType)
54 {
55  // Reasonable to expect any faPatch constraint has an identically
56  // named polyPatch/pointPatch equivalent
57 
58  return polyPatch::constraintType(patchType);
59 }
60 
61 
63 {
64  const auto& cnstrTable = *dictionaryConstructorTablePtr_;
65 
66  wordList cTypes(cnstrTable.size());
67 
68  label i = 0;
69 
70  forAllConstIters(cnstrTable, iter)
71  {
72  if (constraintType(iter.key()))
73  {
74  cTypes[i++] = iter.key();
75  }
76  }
77 
78  cTypes.resize(i);
79 
80  return cTypes;
81 }
82 
83 
84 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
85 
86 void Foam::faPatch::clearOut()
87 {
88  deleteDemandDrivenData(edgeFacesPtr_);
89  deleteDemandDrivenData(pointLabelsPtr_);
90  deleteDemandDrivenData(pointEdgesPtr_);
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
96 Foam::faPatch::faPatch
97 (
98  const word& name,
99  const labelUList& edgeLabels,
100  const label index,
101  const faBoundaryMesh& bm,
102  const label nbrPolyPatchi,
103  const word& patchType
104 )
105 :
106  patchIdentifier(name, index),
107  labelList(edgeLabels),
108  nbrPolyPatchId_(nbrPolyPatchi),
109  boundaryMesh_(bm),
110  edgeFacesPtr_(nullptr),
111  pointLabelsPtr_(nullptr),
112  pointEdgesPtr_(nullptr)
113 {
114  if (constraintType(patchType))
115  {
116  inGroups().appendUniq(patchType);
117  }
118 }
119 
120 
121 Foam::faPatch::faPatch
122 (
123  const word& name,
124  const dictionary& dict,
125  const label index,
126  const faBoundaryMesh& bm,
127  const word& patchType
128 )
129 :
130  patchIdentifier(name, dict, index),
131  labelList(dict.get<labelList>("edgeLabels")),
132  nbrPolyPatchId_(dict.get<label>("ngbPolyPatchIndex")),
133  boundaryMesh_(bm),
134  edgeFacesPtr_(nullptr),
135  pointLabelsPtr_(nullptr),
136  pointEdgesPtr_(nullptr)
137 {
138  if (constraintType(patchType))
139  {
140  inGroups().appendUniq(patchType);
141  }
142 }
143 
144 
145 Foam::faPatch::faPatch
146 (
147  const faPatch& p,
148  const faBoundaryMesh& bm,
149  const label index,
150  const labelUList& edgeLabels,
151  const label nbrPolyPatchi
152 )
153 :
154  patchIdentifier(p, index),
155  labelList(edgeLabels),
156  nbrPolyPatchId_(p.nbrPolyPatchId_),
157  boundaryMesh_(bm),
158  edgeFacesPtr_(nullptr),
159  pointLabelsPtr_(nullptr),
160  pointEdgesPtr_(nullptr)
161 {}
162 
163 
164 Foam::faPatch::faPatch
165 (
166  const faPatch& p,
167  const faBoundaryMesh& bm
168 )
169 :
170  faPatch
171  (
172  p,
173  bm,
174  p.index(),
175  p.edgeLabels(),
176  p.nbrPolyPatchId_
177  )
178 {}
179 
180 
181 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
182 
184 {
185  clearOut();
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 {
193  return boundaryMesh_;
194 }
195 
196 
197 Foam::label Foam::faPatch::offset() const
198 {
199  return max
200  (
201  0,
202  boundaryMesh().mesh().patchStarts()[index()]
203  - boundaryMesh().mesh().nInternalEdges()
204  );
205 }
206 
208 Foam::label Foam::faPatch::start() const
209 {
210  return boundaryMesh().mesh().patchStarts()[index()];
211 }
212 
213 
215 {
216  const auto& connections = boundaryMesh().mesh().boundaryConnections();
217  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
218 
219  List<labelPair> output(this->nEdges());
220 
221  // Like an IndirectList but removing the nInternalEdges offset
222  label count = 0;
223  for (const label patchEdgei : this->edgeLabels())
224  {
225  const label bndEdgei = (patchEdgei - nInternalEdges);
226  output[count] = connections[bndEdgei];
227  ++count;
228  }
229 
230  return output;
231 }
232 
233 
235 {
236  const auto& connections = boundaryMesh().mesh().boundaryConnections();
237  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
238 
239  labelHashSet procsUsed(2*Pstream::nProcs());
240 
241  for (const label patchEdgei : this->edgeLabels())
242  {
243  const label bndEdgei = (patchEdgei - nInternalEdges);
244  const label proci = connections[bndEdgei].first();
245  procsUsed.insert(proci);
246  }
247  procsUsed.erase(-1); // placeholder value
248  procsUsed.erase(Pstream::myProcNo());
249 
250  return procsUsed.sortedToc();
251 }
252 
253 
255 {
256  const auto& connections = boundaryMesh().mesh().boundaryConnections();
257  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
258 
259  Map<label> procCount(2*Pstream::nProcs());
260 
261  for (const label patchEdgei : this->edgeLabels())
262  {
263  const label bndEdgei = (patchEdgei - nInternalEdges);
264  const label proci = connections[bndEdgei].first();
265  ++procCount(proci);
266  }
267  procCount.erase(-1); // placeholder value
268  procCount.erase(Pstream::myProcNo());
269 
270  // Flatten as list
271  List<labelPair> output(procCount.size());
272  label count = 0;
273  for (const label proci : procCount.sortedToc())
274  {
275  output[count].first() = proci;
276  output[count].second() = procCount[proci]; // size
277  ++count;
278  }
279 
280  return output;
281 }
282 
283 
285 {
286  if (!pointLabelsPtr_)
287  {
288  calcPointLabels();
289  }
290 
291  return *pointLabelsPtr_;
292 }
293 
294 
296 {
297  if (!pointEdgesPtr_)
298  {
299  calcPointEdges();
300  }
301 
302  return *pointEdgesPtr_;
303 }
304 
305 
307 {
308  const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
309 
310  // Walk boundary edges.
311  // The edge orientation corresponds to the face orientation
312  // (outwards normal).
313 
314  // Note: could combine this with calcPointEdges for more efficiency
315 
316  // Map<label> markedPoints(4*edges.size());
317  labelHashSet markedPoints(4*edges.size());
318  DynamicList<label> dynEdgePoints(2*edges.size());
319 
320  for (const edge& e : edges)
321  {
322  // if (markedPoints.insert(e.first(), markedPoints.size()))
323  if (markedPoints.insert(e.first()))
324  {
325  dynEdgePoints.append(e.first());
326  }
327  // if (markedPoints.insert(e.second(), markedPoints.size()))
328  if (markedPoints.insert(e.second()))
329  {
330  dynEdgePoints.append(e.second());
331  }
332  }
333 
334  // Transfer to plain list (reuse storage)
335  pointLabelsPtr_ = new labelList(std::move(dynEdgePoints));
360 }
361 
362 
364 {
365  const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
366 
367  const labelList& edgePoints = pointLabels();
368 
369  // Cannot use invertManyToMany - we have non-local edge numbering
370 
371  // Intermediate storage for pointEdges.
372  // Points on the boundary will normally connect 1 or 2 edges only.
373  List<DynamicList<label,2>> dynPointEdges(edgePoints.size());
374 
375  forAll(edges, edgei)
376  {
377  const edge& e = edges[edgei];
378 
379  dynPointEdges[edgePoints.find(e.first())].append(edgei);
380  dynPointEdges[edgePoints.find(e.second())].append(edgei);
381  }
382 
383  // Flatten to regular list
384  pointEdgesPtr_ = new labelListList(edgePoints.size());
385  auto& pEdges = *pointEdgesPtr_;
386 
387  forAll(pEdges, pointi)
388  {
389  pEdges[pointi] = std::move(dynPointEdges[pointi]);
390  }
391 }
392 
393 
395 {
396  if (nbrPolyPatchId_ < 0)
397  {
399  }
400 
401  return boundaryMesh().mesh().haloFaceNormals(this->index());
402 }
403 
404 
406 {
407  if (nbrPolyPatchId_ < 0)
408  {
409  return tmp<vectorField>::New();
410  }
411 
412  // Unit normals for the neighbour patch faces
413  const vectorField faceNormals
414  (
415  boundaryMesh().mesh().haloFaceNormals(this->index())
416  );
417 
418  const labelListList& pntEdges = pointEdges();
419 
420  auto tpointNorm = tmp<vectorField>::New(pntEdges.size());
421  auto& pointNorm = tpointNorm.ref();
422 
423  forAll(pointNorm, pointi)
424  {
425  vector& n = pointNorm[pointi];
426  n = Zero;
427 
428  for (const label bndEdgei : pntEdges[pointi])
429  {
430  n += faceNormals[bndEdgei];
431  }
432 
433  n.normalise();
434  }
435 
436  return tpointNorm;
437 }
438 
439 
441 {
442  if (!edgeFacesPtr_)
443  {
444  edgeFacesPtr_ = new labelList::subList
445  (
446  patchSlice(boundaryMesh().mesh().edgeOwner())
447  );
448  }
449 
450  return *edgeFacesPtr_;
451 }
452 
455 {
456  return boundaryMesh().mesh().edgeCentres().boundaryField()[index()];
457 }
458 
461 {
462  return boundaryMesh().mesh().Le().boundaryField()[index()];
463 }
464 
467 {
468  return boundaryMesh().mesh().magLe().boundaryField()[index()];
469 }
470 
471 
473 {
474  auto tedgeNorm = tmp<vectorField>::New(edgeLengths());
476  tedgeNorm.ref().normalise();
477 
478  return tedgeNorm;
479 }
480 
483 {
484  return patchInternalField(boundaryMesh().mesh().areaCentres());
485 }
486 
487 
489 {
490  // Use patch-normal delta for all non-coupled BCs
491  const vectorField nHat(edgeNormals());
492 
493  vectorField edgePN(edgeCentres() - edgeFaceCentres());
494 
495  // Do not allow any mag(val) < SMALL
496  // sqrt(1/3) = 0.5773502691896257, but slightly rounded down
497  const vector minVector(vector::uniform(0.57735*SMALL));
498  const scalar minLenSqr(SMALL*SMALL);
499 
500  for (vector& e : edgePN)
501  {
502  if (e.magSqr() < minLenSqr)
503  {
504  e = minVector;
505  }
506  }
507 
508  return nHat*(nHat & edgePN);
509 }
510 
513 {
514  dc = scalar(1)/(edgeNormals() & delta());
515 }
516 
517 
519 {
520  vectorField unitDelta(delta()/mag(delta()));
521  vectorField edgeNormMag(edgeNormals()/mag(edgeNormals()));
522  scalarField dn(edgeNormals() & delta());
523 
524  k = edgeNormMag - (scalar(1)/(unitDelta & edgeNormMag))*unitDelta;
525 }
526 
529 {
530  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
531 }
532 
535 {
536  w = scalar(1);
537 }
538 
541 {
542  return boundaryMesh().mesh().weights().boundaryField()[index()];
543 }
544 
545 
547 {}
548 
549 
551 {
552  clearOut();
553  static_cast<labelList&>(*this) = newEdges;
554 }
555 
556 
558 {
559  clearOut();
560  static_cast<labelList&>(*this) = std::move(newEdges);
561 }
562 
563 
564 void Foam::faPatch::write(Ostream& os) const
565 {
566  os.writeEntry("type", type());
567 
569 
570  os.writeEntry("ngbPolyPatchIndex", nbrPolyPatchId_);
571  static_cast<const labelList&>(*this).writeEntry("edgeLabels", os);
572 }
573 
574 
575 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
576 
577 Foam::Ostream& Foam::operator<<(Ostream& os, const faPatch& p)
578 {
579  p.write(os);
581  return os;
582 }
583 
584 
585 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
scalar delta
dictionary dict
const scalarField & weights() const
Return patch weighting factors.
Definition: faPatch.C:533
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:481
void write(Ostream &os) const
Write (physicalType, inGroups) dictionary entries (without surrounding braces)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList pointLabels(nPoints, -1)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
const scalarField & deltaCoeffs() const
Return patch edge - neighbour face distances.
Definition: faPatch.C:521
Identifies a patch by name and index, with optional physical type and group information.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: faPatch.C:46
const wordList & inGroups() const noexcept
The (optional) groups that the patch belongs to.
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:535
const labelList & pointLabels() const
Return patch point labels.
Definition: faPatch.C:277
label nInternalEdges() const
Number of internal edges.
virtual void makeDeltaCoeffs(scalarField &) const
Make patch edge - neighbour face distances.
Definition: faPatch.C:505
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
const scalarField & magEdgeLengths() const
Return edge length magnitudes, like the faMesh::magLe() method.
Definition: faPatch.C:459
SubList< edge > subList
Declare type of subList.
Definition: List.H:122
label k
Boltzmann constant.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1029
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
List< labelPair > boundaryProcSizes() const
List of proc/size for the boundary edge neighbour processors (does not include own proc) ...
Definition: faPatch.C:247
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: faPatch.C:55
void makeCorrectionVectors(vectorField &) const
Definition: faPatch.C:511
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:331
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
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1020
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const vectorField & edgeCentres() const
Return edge centres.
Definition: faPatch.C:447
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
tmp< vectorField > edgeFaceCentres() const
Return neighbour face centres.
Definition: faPatch.C:475
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
void calcPointLabels() const
Calculate patch point labels.
Definition: faPatch.C:299
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual void write(Ostream &) const
Write.
Definition: faPatch.C:557
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:212
Vector< scalar > vector
Definition: vector.H:57
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition: faPatch.C:433
tmp< vectorField > edgeNormals() const
Return edge unit normals, like the faMesh::unitLe() method.
Definition: faPatch.C:465
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
labelList boundaryProcs() const
Boundary edge neighbour processors (does not include own proc)
Definition: faPatch.C:227
defineTypeNameAndDebug(combustionModel, 0)
List< labelPair > boundaryConnections() const
List of proc/face for the boundary edge neighbours in locally reordered edge numbering.
Definition: faPatch.C:207
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:72
label start() const
Patch start in edge list.
Definition: faPatch.C:201
Buffers for inter-processor communications streams (UOPstream, UIPstream).
virtual ~faPatch()
Destructor.
Definition: faPatch.C:176
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:473
tmp< vectorField > ngbPolyPatchPointNormals() const
Return normals of neighbour polyPatch joined points.
Definition: faPatch.C:398
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:76
Template functions to aid in the implementation of demand driven data.
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: faPatch.C:527
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:270
const vectorField & edgeLengths() const
Return edge length vectors, like the faMesh::Le() method.
Definition: faPatch.C:453
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:130
Finite area boundary mesh.
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition: faPatch.C:288
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:58
label n
void writeEntry(Ostream &os) const
Write the UList with its compound type.
Definition: UListIO.C:31
Field< vector > vectorField
Specialisation of Field<T> for vector.
const bMesh & mesh() const
Definition: boundaryMesh.H:259
label offset() const
The offset where this patch starts in the boundary edge list.
Definition: faPatch.C:190
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void deleteDemandDrivenData(DataPtr &dataPtr)
tmp< vectorField > ngbPolyPatchFaceNormals() const
Return normals of neighbour polyPatch faces.
Definition: faPatch.C:387
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
Definition: faPatch.C:539
void calcPointEdges() const
Calculate patch point-edge addressing.
Definition: faPatch.C:356
void resetEdges(const labelUList &newEdges)
Reset the list of edges (use with caution)
Definition: faPatch.C:543
const faBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: faPatch.C:184
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133