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 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(faPatch, 0);
47 }
48 
49 
50 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
51 
52 bool Foam::faPatch::constraintType(const word& patchType)
53 {
54  // Reasonable to expect any faPatch constraint has an identically
55  // named polyPatch/pointPatch equivalent
56 
57  return polyPatch::constraintType(patchType);
58 }
59 
60 
62 {
63  const auto& cnstrTable = *dictionaryConstructorTablePtr_;
64 
65  wordList cTypes(cnstrTable.size());
66 
67  label i = 0;
68 
69  forAllConstIters(cnstrTable, iter)
70  {
71  if (constraintType(iter.key()))
72  {
73  cTypes[i++] = iter.key();
74  }
75  }
76 
77  cTypes.resize(i);
78 
79  return cTypes;
80 }
81 
82 
83 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
84 
85 void Foam::faPatch::clearOut()
86 {
87  edgeFacesPtr_.reset(nullptr);
88  pointLabelsPtr_.reset(nullptr);
89  pointEdgesPtr_.reset(nullptr);
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
95 Foam::faPatch::faPatch
96 (
97  const word& name,
98  const labelUList& edgeLabels,
99  const label index,
100  const faBoundaryMesh& bm,
101  const label nbrPolyPatchi,
102  const word& patchType
103 )
104 :
105  patchIdentifier(name, index),
106  labelList(edgeLabels),
107  nbrPolyPatchId_(nbrPolyPatchi),
108  boundaryMesh_(bm),
109  edgeFacesPtr_(nullptr),
110  pointLabelsPtr_(nullptr),
111  pointEdgesPtr_(nullptr)
112 {
113  if (constraintType(patchType))
114  {
115  addGroup(patchType);
116  }
117 }
118 
119 
120 Foam::faPatch::faPatch
121 (
122  const word& name,
123  const dictionary& dict,
124  const label index,
125  const faBoundaryMesh& bm,
126  const word& patchType
127 )
128 :
129  patchIdentifier(name, dict, index),
130  labelList(dict.get<labelList>("edgeLabels")),
131  nbrPolyPatchId_(dict.get<label>("ngbPolyPatchIndex")),
132  boundaryMesh_(bm),
133  edgeFacesPtr_(nullptr),
134  pointLabelsPtr_(nullptr),
135  pointEdgesPtr_(nullptr)
136 {
137  if (constraintType(patchType))
138  {
139  addGroup(patchType);
140  }
141 }
142 
143 
144 Foam::faPatch::faPatch
145 (
146  const faPatch& p,
147  const faBoundaryMesh& bm,
148  const label index,
149  const labelUList& edgeLabels,
150  const label nbrPolyPatchi
151 )
152 :
153  patchIdentifier(p, index),
154  labelList(edgeLabels),
155  nbrPolyPatchId_(p.nbrPolyPatchId_),
156  boundaryMesh_(bm),
157  edgeFacesPtr_(nullptr),
158  pointLabelsPtr_(nullptr),
159  pointEdgesPtr_(nullptr)
160 {}
161 
162 
163 Foam::faPatch::faPatch
164 (
165  const faPatch& p,
166  const faBoundaryMesh& bm
167 )
168 :
169  faPatch
170  (
171  p,
172  bm,
173  p.index(),
174  p.edgeLabels(),
175  p.nbrPolyPatchId_
176  )
177 {}
178 
179 
180 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
181 
183 {
184  clearOut();
185 }
186 
187 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 {
192  return boundaryMesh_;
193 }
194 
195 
196 Foam::label Foam::faPatch::offset() const
197 {
198  return max
199  (
200  0,
201  boundaryMesh().mesh().patchStarts()[index()]
202  - boundaryMesh().mesh().nInternalEdges()
203  );
204 }
205 
207 Foam::label Foam::faPatch::start() const
208 {
209  return boundaryMesh().mesh().patchStarts()[index()];
210 }
211 
212 
214 {
215  const auto& connections = boundaryMesh().mesh().boundaryConnections();
216  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
217 
218  List<labelPair> output(this->nEdges());
219 
220  // Like an IndirectList but removing the nInternalEdges offset
221  label count = 0;
222  for (const label patchEdgei : this->edgeLabels())
223  {
224  const label bndEdgei = (patchEdgei - nInternalEdges);
225  output[count] = connections[bndEdgei];
226  ++count;
227  }
228 
229  return output;
230 }
231 
232 
234 {
235  const auto& connections = boundaryMesh().mesh().boundaryConnections();
236  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
237 
238  labelHashSet procsUsed(2*Pstream::nProcs());
239 
240  for (const label patchEdgei : this->edgeLabels())
241  {
242  const label bndEdgei = (patchEdgei - nInternalEdges);
243  const label proci = connections[bndEdgei].first();
244  procsUsed.insert(proci);
245  }
246  procsUsed.erase(-1); // placeholder value
247  procsUsed.erase(Pstream::myProcNo());
248 
249  return procsUsed.sortedToc();
250 }
251 
252 
254 {
255  const auto& connections = boundaryMesh().mesh().boundaryConnections();
256  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
257 
258  Map<label> procCount(2*Pstream::nProcs());
259 
260  for (const label patchEdgei : this->edgeLabels())
261  {
262  const label bndEdgei = (patchEdgei - nInternalEdges);
263  const label proci = connections[bndEdgei].first();
264  ++procCount(proci);
265  }
266  procCount.erase(-1); // placeholder value
267  procCount.erase(Pstream::myProcNo());
268 
269  // Flatten as list
270  List<labelPair> output(procCount.size());
271  label count = 0;
272  for (const label proci : procCount.sortedToc())
273  {
274  output[count].first() = proci;
275  output[count].second() = procCount[proci]; // size
276  ++count;
277  }
278 
279  return output;
280 }
281 
282 
284 {
285  if (!pointLabelsPtr_)
286  {
287  calcPointLabels();
288  }
289 
290  return *pointLabelsPtr_;
291 }
292 
293 
295 {
296  if (!pointEdgesPtr_)
297  {
298  calcPointEdges();
299  }
300 
301  return *pointEdgesPtr_;
302 }
303 
304 
306 {
307  const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
308 
309  // Walk boundary edges.
310  // The edge orientation corresponds to the face orientation
311  // (outwards normal).
312 
313  // Note: could combine this with calcPointEdges for more efficiency
314 
315  // Map<label> markedPoints(4*edges.size());
316  labelHashSet markedPoints(4*edges.size());
317  DynamicList<label> dynEdgePoints(2*edges.size());
318 
319  for (const edge& e : edges)
320  {
321  // if (markedPoints.insert(e.first(), markedPoints.size()))
322  if (markedPoints.insert(e.first()))
323  {
324  dynEdgePoints.append(e.first());
325  }
326  // if (markedPoints.insert(e.second(), markedPoints.size()))
327  if (markedPoints.insert(e.second()))
328  {
329  dynEdgePoints.append(e.second());
330  }
331  }
332 
333  // Transfer to plain list (reuse storage)
334  pointLabelsPtr_.reset(new labelList(std::move(dynEdgePoints)));
335 
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_.reset(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_.reset
445  (
447  (
448  patchSlice(boundaryMesh().mesh().edgeOwner())
449  )
450  );
451  }
452 
453  return *edgeFacesPtr_;
454 }
455 
458 {
459  return boundaryMesh().mesh().edgeCentres().boundaryField()[index()];
460 }
461 
464 {
465  return boundaryMesh().mesh().Le().boundaryField()[index()];
466 }
467 
470 {
471  return boundaryMesh().mesh().magLe().boundaryField()[index()];
472 }
473 
474 
476 {
477  auto tedgeNorm = tmp<vectorField>::New(edgeLengths());
479  tedgeNorm.ref().normalise();
480 
481  return tedgeNorm;
482 }
483 
486 {
487  return patchInternalField(boundaryMesh().mesh().areaCentres());
488 }
489 
490 
492 {
493  // Use patch-normal delta for all non-coupled BCs
494  const vectorField nHat(edgeNormals());
495 
496  vectorField edgePN(edgeCentres() - edgeFaceCentres());
497 
498  // Do not allow any mag(val) < SMALL
499  // sqrt(1/3) = 0.5773502691896257, but slightly rounded down
500  const vector minVector(vector::uniform(0.57735*SMALL));
501  const scalar minLenSqr(SMALL*SMALL);
502 
503  for (vector& e : edgePN)
504  {
505  if (e.magSqr() < minLenSqr)
506  {
507  e = minVector;
508  }
509  }
510 
511  return nHat*(nHat & edgePN);
512 }
513 
516 {
517  dc = scalar(1)/(edgeNormals() & delta());
518 }
519 
520 
522 {
523  vectorField unitDelta(delta()/mag(delta()));
524  vectorField edgeNormMag(edgeNormals()/mag(edgeNormals()));
525  scalarField dn(edgeNormals() & delta());
526 
527  k = edgeNormMag - (scalar(1)/(unitDelta & edgeNormMag))*unitDelta;
528 }
529 
532 {
533  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
534 }
535 
538 {
539  w = scalar(1);
540 }
541 
544 {
545  return boundaryMesh().mesh().weights().boundaryField()[index()];
546 }
547 
548 
550 {}
551 
552 
554 {
555  clearOut();
556  static_cast<labelList&>(*this) = newEdges;
557 }
558 
559 
561 {
562  clearOut();
563  static_cast<labelList&>(*this) = std::move(newEdges);
564 }
565 
566 
567 void Foam::faPatch::write(Ostream& os) const
568 {
569  os.writeEntry("type", type());
570 
572 
573  os.writeEntry("ngbPolyPatchIndex", nbrPolyPatchId_);
574  static_cast<const labelList&>(*this).writeEntry("edgeLabels", os);
575 }
576 
577 
578 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
579 
580 Foam::Ostream& Foam::operator<<(Ostream& os, const faPatch& p)
581 {
582  p.write(os);
584  return os;
585 }
586 
587 
588 // ************************************************************************* //
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.
scalar delta
dictionary dict
const scalarField & weights() const
Return patch weighting factors.
Definition: faPatch.C:536
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:484
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:129
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:524
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:45
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
const labelList & pointLabels() const
Return patch point labels.
Definition: faPatch.C:276
label nInternalEdges() const
Number of internal edges.
virtual void makeDeltaCoeffs(scalarField &) const
Make patch edge - neighbour face distances.
Definition: faPatch.C:508
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
const scalarField & magEdgeLengths() const
Return edge length magnitudes, like the faMesh::magLe() method.
Definition: faPatch.C:462
SubList< edge > subList
Declare type of subList.
Definition: List.H:144
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:1074
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:246
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: faPatch.C:54
void makeCorrectionVectors(vectorField &) const
Definition: faPatch.C:514
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:342
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:1065
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const vectorField & edgeCentres() const
Return edge centres.
Definition: faPatch.C:450
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
tmp< vectorField > edgeFaceCentres() const
Return neighbour face centres.
Definition: faPatch.C:478
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:298
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:560
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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:468
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
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:226
defineTypeNameAndDebug(combustionModel, 0)
List< labelPair > boundaryConnections() const
List of proc/face for the boundary edge neighbours in locally reordered edge numbering.
Definition: faPatch.C:206
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:200
Buffers for inter-processor communications streams (UOPstream, UIPstream).
virtual ~faPatch()
Destructor.
Definition: faPatch.C:175
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:472
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:77
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:530
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:269
const vectorField & edgeLengths() const
Return edge length vectors, like the faMesh::Le() method.
Definition: faPatch.C:456
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
Finite area boundary mesh.
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition: faPatch.C:287
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:44
label n
void writeEntry(Ostream &os) const
Write the UList with its compound type.
Definition: UListIO.C:30
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:189
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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:542
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:546
const faBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: faPatch.C:183
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127