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-2021 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 "demandDrivenData.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(faPatch, 0);
44  defineRunTimeSelectionTable(faPatch, dictionary);
45  addToRunTimeSelectionTable(faPatch, faPatch, dictionary);
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::faPatch::clearOut()
52 {
53  deleteDemandDrivenData(edgeFacesPtr_);
54  deleteDemandDrivenData(pointLabelsPtr_);
55  deleteDemandDrivenData(pointEdgesPtr_);
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 Foam::faPatch::faPatch
62 (
63  const word& name,
64  const labelList& edgeLabels,
65  const label index,
66  const faBoundaryMesh& bm,
67  const label ngbPolyPatchIndex
68 )
69 :
70  patchIdentifier(name, index),
71  labelList(edgeLabels),
72  nbrPolyPatchId_(ngbPolyPatchIndex),
73  boundaryMesh_(bm),
74  edgeFacesPtr_(nullptr),
75  pointLabelsPtr_(nullptr),
76  pointEdgesPtr_(nullptr)
77 {}
78 
79 
80 Foam::faPatch::faPatch
81 (
82  const word& name,
83  const dictionary& dict,
84  const label index,
85  const faBoundaryMesh& bm
86 )
87 :
88  patchIdentifier(name, dict, index),
89  labelList(dict.get<labelList>("edgeLabels")),
90  nbrPolyPatchId_(dict.get<label>("ngbPolyPatchIndex")),
91  boundaryMesh_(bm),
92  edgeFacesPtr_(nullptr),
93  pointLabelsPtr_(nullptr),
94  pointEdgesPtr_(nullptr)
95 {}
96 
97 
98 Foam::faPatch::faPatch(const faPatch& p, const faBoundaryMesh& bm)
99 :
100  patchIdentifier(p, p.index()),
101  labelList(p),
102  nbrPolyPatchId_(p.nbrPolyPatchId_),
103  boundaryMesh_(bm),
104  edgeFacesPtr_(nullptr),
105  pointLabelsPtr_(nullptr),
106  pointEdgesPtr_(nullptr)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
111 
113 {
114  clearOut();
115 }
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
121 {
122  return boundaryMesh_;
123 }
124 
125 
126 Foam::label Foam::faPatch::start() const
127 {
128  return boundaryMesh().mesh().patchStarts()[index()];
129 }
130 
131 
133 {
134  const auto& connections = boundaryMesh().mesh().boundaryConnections();
135  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
136 
137  List<labelPair> output(this->nEdges());
138 
139  // Like an IndirectList but removing the nInternalEdges offset
140  label count = 0;
141  for (const label patchEdgei : this->edgeLabels())
142  {
143  const label bndEdgei = (patchEdgei - nInternalEdges);
144  output[count] = connections[bndEdgei];
145  ++count;
146  }
147 
148  return output;
149 }
150 
151 
153 {
154  const auto& connections = boundaryMesh().mesh().boundaryConnections();
155  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
156 
157  labelHashSet procsUsed(2*Pstream::nProcs());
158 
159  for (const label patchEdgei : this->edgeLabels())
160  {
161  const label bndEdgei = (patchEdgei - nInternalEdges);
162  const label proci = connections[bndEdgei].first();
163  procsUsed.insert(proci);
164  }
165  procsUsed.erase(-1); // placeholder value
166  procsUsed.erase(Pstream::myProcNo());
167 
168  return procsUsed.sortedToc();
169 }
170 
171 
173 {
174  const auto& connections = boundaryMesh().mesh().boundaryConnections();
175  const label nInternalEdges = boundaryMesh().mesh().nInternalEdges();
176 
177  Map<label> procCount(2*Pstream::nProcs());
178 
179  for (const label patchEdgei : this->edgeLabels())
180  {
181  const label bndEdgei = (patchEdgei - nInternalEdges);
182  const label proci = connections[bndEdgei].first();
183  ++procCount(proci);
184  }
185  procCount.erase(-1); // placeholder value
186  procCount.erase(Pstream::myProcNo());
187 
188  // Flatten as list
189  List<labelPair> output(procCount.size());
190  label count = 0;
191  for (const label proci : procCount.sortedToc())
192  {
193  output[count].first() = proci;
194  output[count].second() = procCount[proci]; // size
195  ++count;
196  }
197 
198  return output;
199 }
200 
201 
203 {
204  if (!pointLabelsPtr_)
205  {
206  calcPointLabels();
207  }
208 
209  return *pointLabelsPtr_;
210 }
211 
212 
214 {
215  if (!pointEdgesPtr_)
216  {
217  calcPointEdges();
218  }
219 
220  return *pointEdgesPtr_;
221 }
222 
223 
225 {
226  const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
227 
228  // Walk boundary edges.
229  // The edge orientation corresponds to the face orientation
230  // (outwards normal).
231 
232  // Note: could combine this with calcPointEdges for more efficiency
233 
234  // Map<label> markedPoints(4*edges.size());
235  labelHashSet markedPoints(4*edges.size());
236  DynamicList<label> dynEdgePoints(2*edges.size());
237 
238  for (const edge& e : edges)
239  {
240  // if (markedPoints.insert(e.first(), markedPoints.size()))
241  if (markedPoints.insert(e.first()))
242  {
243  dynEdgePoints.append(e.first());
244  }
245  // if (markedPoints.insert(e.second(), markedPoints.size()))
246  if (markedPoints.insert(e.second()))
247  {
248  dynEdgePoints.append(e.second());
249  }
250  }
251 
252  // Transfer to plain list (reuse storage)
253  pointLabelsPtr_ = new labelList(std::move(dynEdgePoints));
278 }
279 
280 
282 {
283  const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
284 
285  const labelList& edgePoints = pointLabels();
286 
287  // Cannot use invertManyToMany - we have non-local edge numbering
288 
289  // Intermediate storage for pointEdges.
290  // Points on the boundary will normally connect 1 or 2 edges only.
291  List<DynamicList<label,2>> dynPointEdges(edgePoints.size());
292 
293  forAll(edges, edgei)
294  {
295  const edge& e = edges[edgei];
296 
297  dynPointEdges[edgePoints.find(e.first())].append(edgei);
298  dynPointEdges[edgePoints.find(e.second())].append(edgei);
299  }
300 
301  // Flatten to regular list
302  pointEdgesPtr_ = new labelListList(edgePoints.size());
303  auto& pEdges = *pointEdgesPtr_;
304 
305  forAll(pEdges, pointi)
306  {
307  pEdges[pointi] = std::move(dynPointEdges[pointi]);
308  }
309 }
310 
311 
313 {
314  if (nbrPolyPatchId_ < 0)
315  {
316  return tmp<vectorField>::New();
317  }
318 
319  return boundaryMesh().mesh().haloFaceNormals(this->index());
320 }
321 
322 
324 {
325  if (nbrPolyPatchId_ < 0)
326  {
327  return tmp<vectorField>::New();
328  }
329 
330  // Unit normals for the neighbour patch faces
332  (
333  boundaryMesh().mesh().haloFaceNormals(this->index())
334  );
335 
336  const labelListList& pntEdges = pointEdges();
337 
338  auto tpointNorm = tmp<vectorField>::New(pntEdges.size());
339  auto& pointNorm = tpointNorm.ref();
340 
341  forAll(pointNorm, pointi)
342  {
343  vector& n = pointNorm[pointi];
344  n = Zero;
345 
346  for (const label bndEdgei : pntEdges[pointi])
347  {
348  n += faceNormals[bndEdgei];
349  }
350 
351  n.normalise();
352  }
353 
354  return tpointNorm;
355 }
356 
357 
359 {
360  if (!edgeFacesPtr_)
361  {
362  edgeFacesPtr_ = new labelList::subList
363  (
364  patchSlice(boundaryMesh().mesh().edgeOwner())
365  );
366  }
367 
368  return *edgeFacesPtr_;
369 }
370 
371 
373 {
374  return boundaryMesh().mesh().edgeCentres().boundaryField()[index()];
375 }
376 
377 
379 {
380  return boundaryMesh().mesh().Le().boundaryField()[index()];
381 }
382 
383 
385 {
386  return boundaryMesh().mesh().magLe().boundaryField()[index()];
387 }
388 
389 
391 {
392  auto tedgeNorm = tmp<vectorField>::New(edgeLengths());
393 
394  for (vector& n : tedgeNorm.ref())
395  {
396  n.normalise();
397  }
398 
399  return tedgeNorm;
400 }
401 
402 
404 {
405  auto tfc = tmp<vectorField>::New(size());
406  auto& fc = tfc.ref();
407 
408  // get reference to global face centres
409  const vectorField& gfc =
410  boundaryMesh().mesh().areaCentres().internalField();
411 
412  const labelUList& faceLabels = edgeFaces();
413 
414  forAll(faceLabels, edgeI)
415  {
416  fc[edgeI] = gfc[faceLabels[edgeI]];
417  }
418 
419  return tfc;
420 }
421 
422 
424 {
425  return edgeNormals()*(edgeNormals() & (edgeCentres() - edgeFaceCentres()));
426 }
427 
428 
430 {
431  dc = scalar(1)/(edgeNormals() & delta());
432 }
433 
434 
436 {
437  vectorField unitDelta(delta()/mag(delta()));
438  vectorField edgeNormMag(edgeNormals()/mag(edgeNormals()));
439  scalarField dn(edgeNormals() & delta());
440 
441  k = edgeNormMag - (scalar(1)/(unitDelta & edgeNormMag))*unitDelta;
442 }
443 
444 
446 {
447  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
448 }
449 
450 
452 {
453  w = scalar(1);
454 }
455 
456 
458 {
459  return boundaryMesh().mesh().weights().boundaryField()[index()];
460 }
461 
462 
464 {}
465 
466 
468 {
469  clearOut();
470  static_cast<labelList&>(*this) = newEdges;
471 }
472 
473 
475 {
476  clearOut();
477  static_cast<labelList&>(*this) = std::move(newEdges);
478 }
479 
480 
482 {
483  os.writeEntry("type", type());
484 
486 
487  os.writeEntry("ngbPolyPatchIndex", nbrPolyPatchId_);
488  static_cast<const labelList&>(*this).writeEntry("edgeLabels", os);
489 }
490 
491 
492 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
493 
495 {
496  p.write(os);
498  return os;
499 }
500 
501 
502 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:206
Foam::faPatch::edgeFaceCentres
tmp< vectorField > edgeFaceCentres() const
Return neighbour face centres.
Definition: faPatch.C:403
Foam::faPatch::pointEdges
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition: faPatch.C:213
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::output
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::faPatch::movePoints
virtual void movePoints(const pointField &)
Correct patch after moving points.
Definition: faPatch.C:463
Foam::faPatch::resetEdges
void resetEdges(const UList< label > &newEdges)
Reset the list of edges (use with caution)
Definition: faPatch.C:467
Foam::DynamicList< label >
Foam::faPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Return patch edge - neighbour face distances.
Definition: faPatch.C:445
Foam::faPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: faPatch.C:451
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::Map< label >
Foam::faPatch::calcPointLabels
void calcPointLabels() const
Calculate patch point labels.
Definition: faPatch.C:224
faMesh.H
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
append
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::faPatch::boundaryMesh
const faBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: faPatch.C:120
Foam::faPatch::delta
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:423
Foam::faPatch::weights
const scalarField & weights() const
Return patch weighting factors.
Definition: faPatch.C:457
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::faPatch::ngbPolyPatchFaceNormals
tmp< vectorField > ngbPolyPatchFaceNormals() const
Return normals of neighbour polyPatch faces.
Definition: faPatch.C:312
Foam::faPatch::edgeFaces
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition: faPatch.C:358
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::faBoundaryMesh
Finite area boundary mesh.
Definition: faBoundaryMesh.H:65
Foam::patchIdentifier
Identifies a patch by name and index, with optional physical type and group information.
Definition: patchIdentifier.H:55
Foam::faPatch::boundaryProcs
labelList boundaryProcs() const
Definition: faPatch.C:152
faceNormals
surfaceVectorField faceNormals(mesh.Sf()/mesh.magSf())
Foam::Field< vector >
Foam::faPatch::makeDeltaCoeffs
virtual void makeDeltaCoeffs(scalarField &) const
Make patch edge - neighbour face distances.
Definition: faPatch.C:429
edgeFields.H
Foam::faPatch::edgeCentres
const vectorField & edgeCentres() const
Return edge centres.
Definition: faPatch.C:372
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
areaFields.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::List::subList
SubList< T > subList
Declare type of subList.
Definition: List.H:112
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
faBoundaryMesh.H
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::faPatch::boundaryProcSizes
List< labelPair > boundaryProcSizes() const
Definition: faPatch.C:172
Foam::faPatch::makeCorrectionVectors
void makeCorrectionVectors(vectorField &) const
Definition: faPatch.C:435
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faPatch::start
label start() const
Patch start in edge list.
Definition: faPatch.C:126
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatch.C:214
Foam::faPatch::~faPatch
virtual ~faPatch()
Destructor.
Definition: faPatch.C:112
Foam::faPatch::boundaryConnections
List< labelPair > boundaryConnections() const
Definition: faPatch.C:132
Foam::faPatch::pointLabels
const labelList & pointLabels() const
Return patch point labels.
Definition: faPatch.C:202
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
edgeHashes.H
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::faPatch::magEdgeLengths
const scalarField & magEdgeLengths() const
Return edge length magnitudes.
Definition: faPatch.C:384
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::faPatch::write
virtual void write(Ostream &) const
Write.
Definition: faPatch.C:481
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
faPatch.H
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::faPatch::edgeNormals
tmp< vectorField > edgeNormals() const
Return edge normals.
Definition: faPatch.C:390
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::patchIdentifier::write
void write(Ostream &os) const
Definition: patchIdentifier.C:96
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::faPatch::ngbPolyPatchPointNormals
tmp< vectorField > ngbPolyPatchPointNormals() const
Return normals of neighbour polyPatch joined points.
Definition: faPatch.C:323
Foam::faPatch
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:69
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointLabels
labelList pointLabels(nPoints, -1)
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::faPatch::edgeLengths
const vectorField & edgeLengths() const
Return edge length vectors.
Definition: faPatch.C:378
Foam::faPatch::calcPointEdges
void calcPointEdges() const
Calculate patch point-edge addressing.
Definition: faPatch.C:281