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 {
110  if (constraintType(patchType))
111  {
112  addGroup(patchType);
113  }
114 }
115 
116 
117 Foam::faPatch::faPatch
118 (
119  const word& name,
120  const dictionary& dict,
121  const label index,
122  const faBoundaryMesh& bm,
123  const word& patchType
124 )
125 :
126  patchIdentifier(name, dict, index),
127  labelList(dict.get<labelList>("edgeLabels")),
128  nbrPolyPatchId_(dict.get<label>("ngbPolyPatchIndex")),
129  boundaryMesh_(bm)
130 {
131  if (constraintType(patchType))
132  {
133  addGroup(patchType);
134  }
135 }
136 
137 
138 Foam::faPatch::faPatch
139 (
140  const faPatch& p,
141  const faBoundaryMesh& bm,
142  const label index,
143  const labelUList& edgeLabels,
144  const label nbrPolyPatchi
145 )
146 :
147  patchIdentifier(p, index),
148  labelList(edgeLabels),
149  nbrPolyPatchId_(p.nbrPolyPatchId_),
150  boundaryMesh_(bm)
151 {}
152 
153 
154 Foam::faPatch::faPatch
155 (
156  const faPatch& p,
157  const faBoundaryMesh& bm
158 )
159 :
160  faPatch
161  (
162  p,
163  bm,
164  p.index(),
165  p.edgeLabels(),
166  p.nbrPolyPatchId_
167  )
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
172 
174 {
175  clearOut();
176 }
177 
178 
179 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 {
183  return boundaryMesh_;
184 }
185 
186 
187 Foam::label Foam::faPatch::offset() const
188 {
189  return max
190  (
191  0,
192  boundaryMesh().mesh().patchStarts()[index()]
193  - boundaryMesh().mesh().nInternalEdges()
194  );
195 }
196 
198 Foam::label Foam::faPatch::start() const
199 {
200  return boundaryMesh().mesh().patchStarts()[index()];
201 }
202 
203 
205 {
206  const auto& connections = boundaryMesh().mesh().boundaryConnections();
207  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
208 
209  List<labelPair> output(this->nEdges());
210 
211  // Like an IndirectList but removing the nInternalEdges offset
212  label count = 0;
213  for (const label patchEdgei : this->edgeLabels())
214  {
215  const label bndEdgei = (patchEdgei - nInternalEdges);
216  output[count] = connections[bndEdgei];
217  ++count;
218  }
219 
220  return output;
221 }
222 
223 
225 {
226  const auto& connections = boundaryMesh().mesh().boundaryConnections();
227  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
228 
229  labelHashSet procsUsed(2*Pstream::nProcs());
230 
231  for (const label patchEdgei : this->edgeLabels())
232  {
233  const label bndEdgei = (patchEdgei - nInternalEdges);
234  const label proci = connections[bndEdgei].first();
235  procsUsed.insert(proci);
236  }
237  procsUsed.erase(-1); // placeholder value
238  procsUsed.erase(Pstream::myProcNo());
239 
240  return procsUsed.sortedToc();
241 }
242 
243 
245 {
246  const auto& connections = boundaryMesh().mesh().boundaryConnections();
247  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
248 
249  Map<label> procCount(2*Pstream::nProcs());
250 
251  for (const label patchEdgei : this->edgeLabels())
252  {
253  const label bndEdgei = (patchEdgei - nInternalEdges);
254  const label proci = connections[bndEdgei].first();
255  ++procCount(proci);
256  }
257  procCount.erase(-1); // placeholder value
258  procCount.erase(Pstream::myProcNo());
259 
260  // Flatten as list
261  List<labelPair> output(procCount.size());
262  label count = 0;
263  for (const label proci : procCount.sortedToc())
264  {
265  output[count].first() = proci;
266  output[count].second() = procCount[proci]; // size
267  ++count;
268  }
269 
270  return output;
271 }
272 
273 
275 {
276  if (!pointLabelsPtr_)
277  {
278  calcPointLabels();
279  }
280 
281  return *pointLabelsPtr_;
282 }
283 
284 
286 {
287  if (!pointEdgesPtr_)
288  {
289  calcPointEdges();
290  }
291 
292  return *pointEdgesPtr_;
293 }
294 
295 
297 {
298  const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
299 
300  // Walk boundary edges.
301  // The edge orientation corresponds to the face orientation
302  // (outwards normal).
303 
304  // Note: could combine this with calcPointEdges for more efficiency
305 
306  // Map<label> markedPoints(4*edges.size());
307  labelHashSet markedPoints(4*edges.size());
308  DynamicList<label> dynEdgePoints(2*edges.size());
309 
310  for (const edge& e : edges)
311  {
312  // if (markedPoints.insert(e.first(), markedPoints.size()))
313  if (markedPoints.insert(e.first()))
314  {
315  dynEdgePoints.append(e.first());
316  }
317  // if (markedPoints.insert(e.second(), markedPoints.size()))
318  if (markedPoints.insert(e.second()))
319  {
320  dynEdgePoints.append(e.second());
321  }
322  }
323 
324  // Transfer to plain list (reuse storage)
325  pointLabelsPtr_ = std::make_unique<labelList>(std::move(dynEdgePoints));
326 
351 }
352 
353 
355 {
356  const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
357 
358  const labelList& edgePoints = pointLabels();
359 
360  // Cannot use invertManyToMany - we have non-local edge numbering
361 
362  // Intermediate storage for pointEdges.
363  // Points on the boundary will normally connect 1 or 2 edges only.
364  List<DynamicList<label,2>> dynPointEdges(edgePoints.size());
365 
366  forAll(edges, edgei)
367  {
368  const edge& e = edges[edgei];
369 
370  dynPointEdges[edgePoints.find(e.first())].append(edgei);
371  dynPointEdges[edgePoints.find(e.second())].append(edgei);
372  }
373 
374  // Flatten to regular list
375  pointEdgesPtr_ = std::make_unique<labelListList>(edgePoints.size());
376  auto& pEdges = *pointEdgesPtr_;
377 
378  forAll(pEdges, pointi)
379  {
380  pEdges[pointi] = std::move(dynPointEdges[pointi]);
381  }
382 }
383 
384 
386 {
387  if (nbrPolyPatchId_ < 0)
388  {
390  }
391 
392  return boundaryMesh().mesh().haloFaceNormals(this->index());
393 }
394 
395 
397 {
398  if (nbrPolyPatchId_ < 0)
399  {
400  return tmp<vectorField>::New();
401  }
402 
403  // Unit normals for the neighbour patch faces
404  const vectorField faceNormals
405  (
406  boundaryMesh().mesh().haloFaceNormals(this->index())
407  );
408 
409  const labelListList& pntEdges = pointEdges();
410 
411  auto tpointNorm = tmp<vectorField>::New(pntEdges.size());
412  auto& pointNorm = tpointNorm.ref();
413 
414  forAll(pointNorm, pointi)
415  {
416  vector& n = pointNorm[pointi];
417  n = Zero;
418 
419  for (const label bndEdgei : pntEdges[pointi])
420  {
421  n += faceNormals[bndEdgei];
422  }
423 
424  n.normalise();
425  }
426 
427  return tpointNorm;
428 }
429 
430 
432 {
433  if (!edgeFacesPtr_)
434  {
435  edgeFacesPtr_ = std::make_unique<labelList::subList>
436  (
437  patchSlice(boundaryMesh().mesh().edgeOwner())
438  );
439  }
440 
441  return *edgeFacesPtr_;
442 }
443 
446 {
447  return boundaryMesh().mesh().edgeCentres().boundaryField()[index()];
448 }
449 
452 {
453  return boundaryMesh().mesh().Le().boundaryField()[index()];
454 }
455 
458 {
459  return boundaryMesh().mesh().magLe().boundaryField()[index()];
460 }
461 
462 
464 {
465  auto tedgeNorm = tmp<vectorField>::New(edgeLengths());
467  tedgeNorm.ref().normalise();
468 
469  return tedgeNorm;
470 }
471 
474 {
475  return patchInternalField(boundaryMesh().mesh().areaCentres());
476 }
477 
478 
480 {
481  // Use patch-normal delta for all non-coupled BCs
482  const vectorField nHat(edgeNormals());
483 
484  vectorField edgePN(edgeCentres() - edgeFaceCentres());
485 
486  // Do not allow any mag(val) < SMALL
487  // sqrt(1/3) = 0.5773502691896257, but slightly rounded down
488  const vector minVector(vector::uniform(0.57735*SMALL));
489  const scalar minLenSqr(SMALL*SMALL);
490 
491  for (vector& e : edgePN)
492  {
493  if (e.magSqr() < minLenSqr)
494  {
495  e = minVector;
496  }
497  }
498 
499  return nHat*(nHat & edgePN);
500 }
501 
504 {
505  dc = scalar(1)/(edgeNormals() & delta());
506 }
507 
508 
510 {
511  const vectorField unitDelta(delta()/mag(delta()));
512 
513  k = edgeNormals() - (scalar(1)/(unitDelta & edgeNormals()))*unitDelta;
514 }
515 
518 {
519  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
520 }
521 
524 {
525  w = scalar(1);
526 }
527 
530 {
531  return boundaryMesh().mesh().weights().boundaryField()[index()];
532 }
533 
534 
536 {}
537 
538 
540 {
541  clearOut();
542  static_cast<labelList&>(*this) = newEdges;
543 }
544 
545 
547 {
548  clearOut();
549  static_cast<labelList&>(*this) = std::move(newEdges);
550 }
551 
552 
553 void Foam::faPatch::write(Ostream& os) const
554 {
555  os.writeEntry("type", type());
556 
558 
559  os.writeEntry("ngbPolyPatchIndex", nbrPolyPatchId_);
560  static_cast<const labelList&>(*this).writeEntry("edgeLabels", os);
561 }
562 
563 
564 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
565 
566 Foam::Ostream& Foam::operator<<(Ostream& os, const faPatch& p)
567 {
568  p.write(os);
570  return os;
571 }
572 
573 
574 // ************************************************************************* //
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:522
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:472
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:510
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:267
label nInternalEdges() const
Number of internal edges.
virtual void makeDeltaCoeffs(scalarField &) const
Make patch edge - neighbour face distances.
Definition: faPatch.C:496
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:450
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:1086
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:237
#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:502
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:358
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:1077
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const vectorField & edgeCentres() const
Return edge centres.
Definition: faPatch.C:438
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:466
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:289
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:546
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:424
tmp< vectorField > edgeNormals() const
Return edge unit normals, like the faMesh::unitLe() method.
Definition: faPatch.C:456
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:217
defineTypeNameAndDebug(combustionModel, 0)
List< labelPair > boundaryConnections() const
List of proc/face for the boundary edge neighbours in locally reordered edge numbering.
Definition: faPatch.C:197
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:191
Buffers for inter-processor communications streams (UOPstream, UIPstream).
virtual ~faPatch()
Destructor.
Definition: faPatch.C:166
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:496
tmp< vectorField > ngbPolyPatchPointNormals() const
Return normals of neighbour polyPatch joined points.
Definition: faPatch.C:389
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:516
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:255
const vectorField & edgeLengths() const
Return edge length vectors, like the faMesh::Le() method.
Definition: faPatch.C:444
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:163
Finite area boundary mesh.
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition: faPatch.C:278
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:180
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:378
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
Definition: faPatch.C:528
void calcPointEdges() const
Calculate patch point-edge addressing.
Definition: faPatch.C:347
void resetEdges(const labelUList &newEdges)
Reset the list of edges (use with caution)
Definition: faPatch.C:532
const faBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: faPatch.C:174
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