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 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 "polyMesh.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(faPatch, 0);
42  defineRunTimeSelectionTable(faPatch, dictionary);
43  addToRunTimeSelectionTable(faPatch, faPatch, dictionary);
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::faPatch::clearOut()
50 {
51  deleteDemandDrivenData(edgeFacesPtr_);
52  deleteDemandDrivenData(pointLabelsPtr_);
53  deleteDemandDrivenData(pointEdgesPtr_);
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
59 Foam::faPatch::faPatch
60 (
61  const word& name,
62  const labelList& edgeLabels,
63  const label index,
64  const faBoundaryMesh& bm,
65  const label ngbPolyPatchIndex
66 )
67 :
68  labelList(edgeLabels),
69  patchIdentifier(name, index),
70  ngbPolyPatchIndex_(ngbPolyPatchIndex),
71  boundaryMesh_(bm),
72  edgeFacesPtr_(nullptr),
73  pointLabelsPtr_(nullptr),
74  pointEdgesPtr_(nullptr)
75 {}
76 
77 
78 Foam::faPatch::faPatch
79 (
80  const word& name,
81  const dictionary& dict,
82  const label index,
83  const faBoundaryMesh& bm
84 )
85 :
86  labelList(dict.lookup("edgeLabels")),
87  patchIdentifier(name, dict, index),
88  ngbPolyPatchIndex_(dict.get<label>("ngbPolyPatchIndex")),
89  boundaryMesh_(bm),
90  edgeFacesPtr_(nullptr),
91  pointLabelsPtr_(nullptr),
92  pointEdgesPtr_(nullptr)
93 {}
94 
95 
96 Foam::faPatch::faPatch(const faPatch& p, const faBoundaryMesh& bm)
97 :
98  labelList(p),
99  patchIdentifier(p, p.index()),
100  ngbPolyPatchIndex_(p.ngbPolyPatchIndex_),
101  boundaryMesh_(bm),
102  edgeFacesPtr_(nullptr),
103  pointLabelsPtr_(nullptr),
104  pointEdgesPtr_(nullptr)
105 {}
106 
107 
108 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
109 
111 {
112  clearOut();
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  return ngbPolyPatchIndex_;
121 }
122 
123 
125 {
126  return boundaryMesh_;
127 }
128 
129 
131 {
132  return boundaryMesh().mesh().patchStarts()[index()];
133 }
134 
135 
137 {
138  if (!pointLabelsPtr_)
139  {
140  calcPointLabels();
141  }
142 
143  return *pointLabelsPtr_;
144 }
145 
146 
148 {
149  SLList<label> labels;
150 
151  UList<edge> edges = patchSlice(boundaryMesh().mesh().edges());
152 
153  forAll(edges, edgeI)
154  {
155  bool existStart = false;
156  bool existEnd = false;
157 
158  forAllIters(labels, iter)
159  {
160  if (*iter == edges[edgeI].start())
161  {
162  existStart = true;
163  }
164 
165  if (*iter == edges[edgeI].end())
166  {
167  existEnd = true;
168  }
169  }
170 
171  if (!existStart)
172  {
173  labels.append(edges[edgeI].start());
174  }
175 
176  if (!existEnd)
177  {
178  labels.append(edges[edgeI].end());
179  }
180  }
181 
182  pointLabelsPtr_ = new labelList(labels);
183 }
184 
185 
187 {
188  const labelList& points = pointLabels();
189 
190  const edgeList::subList e = patchSlice(boundaryMesh().mesh().edges());
191 
192  // set up storage for pointEdges
193  List<SLList<label>> pointEdgs(points.size());
194 
195  forAll(e, edgeI)
196  {
197  const edge& curPoints = e[edgeI];
198 
199  forAll(curPoints, pointI)
200  {
201  const label localPointIndex = points.find(curPoints[pointI]);
202 
203  pointEdgs[localPointIndex].append(edgeI);
204  }
205  }
206 
207  // sort out the list
208  pointEdgesPtr_ = new labelListList(pointEdgs.size());
209  labelListList& pEdges = *pointEdgesPtr_;
210 
211  forAll(pointEdgs, pointI)
212  {
213  pEdges[pointI].setSize(pointEdgs[pointI].size());
214 
215  label i = 0;
216  for
217  (
218  SLList<label>::iterator curEdgesIter = pointEdgs[pointI].begin();
219  curEdgesIter != pointEdgs[pointI].end();
220  ++curEdgesIter, ++i
221  )
222  {
223  pEdges[pointI][i] = curEdgesIter();
224  }
225  }
226 }
227 
228 
230 {
231  if (!pointEdgesPtr_)
232  {
233  calcPointEdges();
234  }
235 
236  return *pointEdgesPtr_;
237 }
238 
239 
241 {
242  labelList ngbFaces;
243 
244  if (ngbPolyPatchIndex() == -1)
245  {
246  return ngbFaces;
247  }
248 
249  ngbFaces.setSize(faPatch::size());
250 
251  const faMesh& aMesh = boundaryMesh().mesh();
252  const polyMesh& pMesh = aMesh();
254 
255  const labelListList& edgeFaces = pMesh.edgeFaces();
256 
257  labelList faceCells(patch.size(), -1);
258 
259  forAll(faceCells, faceI)
260  {
261  label faceID = aMesh.faceLabels()[faceI];
262 
263  faceCells[faceI] = pMesh.faceOwner()[faceID];
264  }
265 
266  labelList meshEdges
267  (
268  patch.meshEdges
269  (
270  pMesh.edges(),
271  pMesh.cellEdges(),
272  faceCells
273  )
274  );
275 
276  forAll(ngbFaces, edgeI)
277  {
278  ngbFaces[edgeI] = -1;
279 
280  label curEdge = (*this)[edgeI];
281 
282  label curPMeshEdge = meshEdges[curEdge];
283 
284  forAll(edgeFaces[curPMeshEdge], faceI)
285  {
286  label curFace = edgeFaces[curPMeshEdge][faceI];
287 
288  label curPatchID = pMesh.boundaryMesh().whichPatch(curFace);
289 
290  if (curPatchID == ngbPolyPatchIndex())
291  {
292  ngbFaces[edgeI] = curFace;
293  }
294  }
295 
296  if (ngbFaces[edgeI] == -1)
297  {
299  << "Problem with determination of edge ngb faces!"
300  << endl;
301  }
302  }
303 
304  return ngbFaces;
305 }
306 
307 
309 {
310  auto tfN = tmp<vectorField>::New();
311  auto& fN = tfN.ref();
312 
313  if (ngbPolyPatchIndex() == -1)
314  {
315  return tfN;
316  }
317 
318  fN.setSize(faPatch::size());
319 
320  labelList ngbFaces = ngbPolyPatchFaces();
321 
322  const polyMesh& pMesh = boundaryMesh().mesh()();
323 
324  const faceList& faces = pMesh.faces();
325  const pointField& points = pMesh.points();
326 
327  forAll(fN, faceI)
328  {
329  fN[faceI] = faces[ngbFaces[faceI]].unitNormal(points);
330  }
331 
332  return tfN;
333 }
334 
335 
337 {
338  if (ngbPolyPatchIndex() == -1)
339  {
340  return tmp<vectorField>::New();
341  }
342 
343  const labelListList& pntEdges = pointEdges();
344 
345  tmp<vectorField> tpN(new vectorField(pntEdges.size(), Zero));
346  vectorField& pN = tpN.ref();
347 
348  const vectorField faceNormals(ngbPolyPatchFaceNormals());
349 
350  forAll(pN, pointI)
351  {
352  forAll(pntEdges[pointI], edgeI)
353  {
354  pN[pointI] += faceNormals[pntEdges[pointI][edgeI]];
355  }
356  }
357 
358  pN /= mag(pN);
359 
360  return tpN;
361 }
362 
363 
365 {
366  if (!edgeFacesPtr_)
367  {
368  edgeFacesPtr_ = new labelList::subList
369  (
370  patchSlice(boundaryMesh().mesh().edgeOwner())
371  );
372  }
373 
374  return *edgeFacesPtr_;
375 }
376 
377 
379 {
380  return boundaryMesh().mesh().edgeCentres().boundaryField()[index()];
381 }
382 
383 
385 {
386  return boundaryMesh().mesh().Le().boundaryField()[index()];
387 }
388 
389 
391 {
392  return boundaryMesh().mesh().magLe().boundaryField()[index()];
393 }
394 
395 
397 {
398  tmp<vectorField> eN(new vectorField(size()));
399 
400  eN.ref() = edgeLengths()/magEdgeLengths();
401 
402  return eN;
403 }
404 
405 
407 {
408  tmp<vectorField> tfc(new vectorField(size()));
409  vectorField& fc = tfc.ref();
410 
411  // get reference to global face centres
412  const vectorField& gfc =
413  boundaryMesh().mesh().areaCentres().internalField();
414 
415  const labelUList& faceLabels = edgeFaces();
416 
417  forAll(faceLabels, edgeI)
418  {
419  fc[edgeI] = gfc[faceLabels[edgeI]];
420  }
421 
422  return tfc;
423 }
424 
425 
427 {
428  return edgeNormals()*(edgeNormals() & (edgeCentres() - edgeFaceCentres()));
429  //return edgeCentres() - edgeFaceCentres();
430 }
431 
432 
434 {
435  dc = 1.0/(edgeNormals() & delta());
436 }
437 
438 
440 {
441  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
442 }
443 
444 
446 {
447  w = 1.0;
448 }
449 
450 
452 {
453  return boundaryMesh().mesh().weights().boundaryField()[index()];
454 }
455 
456 
458 {}
459 
460 
461 void Foam::faPatch::resetEdges(const labelList& newEdges)
462 {
463  Info<< "Resetting patch edges" << endl;
464  labelList::operator=(newEdges);
465 
466  clearOut();
467 }
468 
469 
470 void Foam::faPatch::write(Ostream& os) const
471 {
472  os.writeEntry("type", type());
473 
475 
476  const labelList& edgeLabels = *this;
477  edgeLabels.writeEntry("edgeLabels", os);
478  os.writeEntry("ngbPolyPatchIndex", ngbPolyPatchIndex_);
479 }
480 
481 
482 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
483 
485 {
486  p.write(os);
487  os.check(FUNCTION_NAME);
488  return os;
489 }
490 
491 
492 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::LList::append
void append(const T &item)
Add copy at tail of list.
Definition: LList.H:237
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:205
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::faPatch::ngbPolyPatchIndex
label ngbPolyPatchIndex() const
Return neighbour polyPatch index.
Definition: faPatch.C:119
Foam::faPatch::edgeFaceCentres
tmp< vectorField > edgeFaceCentres() const
Return neighbour face centres.
Definition: faPatch.C:407
Foam::faPatch::pointEdges
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition: faPatch.C:230
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
stdFoam::begin
constexpr auto begin(C &c) -> decltype(c.begin())
Return iterator to the beginning of the container c.
Definition: stdFoam.H:91
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::faMesh::patch
const indirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMesh.C:911
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::faPatch::movePoints
virtual void movePoints(const pointField &)
Correct patch after moving points.
Definition: faPatch.C:458
Foam::faPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Return patch edge - neighbour face distances.
Definition: faPatch.C:440
Foam::faPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: faPatch.C:446
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::primitiveMesh::edgeFaces
const labelListList & edgeFaces() const
Definition: primitiveMeshEdgeFaces.C:33
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
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::faMesh::faceLabels
const labelList & faceLabels() const
Return faMesh face labels.
Definition: faMesh.H:424
Foam::faPatch::calcPointLabels
void calcPointLabels() const
Calculate patch point labels.
Definition: faPatch.C:148
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
faMesh.H
Foam::LList::iterator
An STL-conforming iterator.
Definition: LList.H:322
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
polyMesh.H
Foam::faPatch::delta
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:427
Foam::faPatch::weights
const scalarField & weights() const
Return patch weighting factors.
Definition: faPatch.C:452
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::faPatch::ngbPolyPatchFaceNormals
tmp< vectorField > ngbPolyPatchFaceNormals() const
Return normals of neighbour polyPatch faces.
Definition: faPatch.C:309
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
aMesh
faMesh aMesh(mesh)
Foam::faPatch::edgeFaces
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition: faPatch.C:365
Foam::LList
Template class for non-intrusive linked lists.
Definition: LList.H:54
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::faBoundaryMesh
Finite area boundary mesh.
Definition: faBoundaryMesh.H:66
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::faPatch::boundaryMesh
const faBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: faPatch.C:125
Foam::patchIdentifier
Identifies a patch by name, patch index and physical type.
Definition: patchIdentifier.H:57
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
faceNormals
surfaceVectorField faceNormals(mesh.Sf()/mesh.magSf())
Foam::Field< vector >
Foam::primitiveMesh::cellEdges
const labelListList & cellEdges() const
Definition: primitiveMeshCellEdges.C:120
Foam::faPatch::makeDeltaCoeffs
virtual void makeDeltaCoeffs(scalarField &) const
Make patch edge - neighbour face distances.
Definition: faPatch.C:434
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
edgeFields.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::faPatch::edgeCentres
const vectorField & edgeCentres() const
Return edge centres.
Definition: faPatch.C:379
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:805
Foam::faPatch::ngbPolyPatchFaces
labelList ngbPolyPatchFaces() const
Return edge neighbour polyPatch faces.
Definition: faPatch.C:241
Foam::List< label >::operator=
void operator=(const UList< label > &a)
Assignment to UList operator. Takes linear time.
Definition: List.C:478
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
areaFields.H
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:217
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:419
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::List< label >::subList
SubList< label > subList
Declare type of subList.
Definition: List.H:117
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
faBoundaryMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:115
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faPatch::start
label start() const
Patch start in edge list.
Definition: faPatch.C:131
Foam::faPatch::~faPatch
virtual ~faPatch()
Destructor.
Definition: faPatch.C:111
Foam::faPatch::pointLabels
const labelList & pointLabels() const
Return patch point labels.
Definition: faPatch.C:137
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
Foam::faPatch::magEdgeLengths
const scalarField & magEdgeLengths() const
Return edge length magnitudes.
Definition: faPatch.C:391
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::faPatch::write
virtual void write(Ostream &) const
Write.
Definition: faPatch.C:471
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
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
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:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
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::faPatch::edgeNormals
tmp< vectorField > edgeNormals() const
Return edge normals.
Definition: faPatch.C:397
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:261
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:77
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
Foam::patchIdentifier::write
void write(Ostream &os) const
Write patchIdentifier as a dictionary.
Definition: patchIdentifier.C:85
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::faPatch::ngbPolyPatchPointNormals
tmp< vectorField > ngbPolyPatchPointNormals() const
Return normals of neighbour polyPatch joined points.
Definition: faPatch.C:337
Foam::faPatch
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:65
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::faPatch::resetEdges
void resetEdges(const labelList &)
Reset edge list.
Definition: faPatch.C:462
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
pointLabels
labelList pointLabels(nPoints, -1)
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &)
Definition: boundaryPatch.C:102
Foam::faPatch::edgeLengths
const vectorField & edgeLengths() const
Return edge length vectors.
Definition: faPatch.C:385
Foam::faPatch::calcPointEdges
void calcPointEdges() const
Calculate patch point-edge addressing.
Definition: faPatch.C:187
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90
Foam::faPatch::size
virtual label size() const
Patch size.
Definition: faPatch.H:247