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-2020 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 #include "demandDrivenData.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(faPatch, 0);
43  defineRunTimeSelectionTable(faPatch, dictionary);
44  addToRunTimeSelectionTable(faPatch, faPatch, dictionary);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::faPatch::clearOut()
51 {
52  deleteDemandDrivenData(edgeFacesPtr_);
53  deleteDemandDrivenData(pointLabelsPtr_);
54  deleteDemandDrivenData(pointEdgesPtr_);
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
60 Foam::faPatch::faPatch
61 (
62  const word& name,
63  const labelList& edgeLabels,
64  const label index,
65  const faBoundaryMesh& bm,
66  const label ngbPolyPatchIndex
67 )
68 :
69  labelList(edgeLabels),
70  patchIdentifier(name, index),
71  ngbPolyPatchIndex_(ngbPolyPatchIndex),
72  boundaryMesh_(bm),
73  edgeFacesPtr_(nullptr),
74  pointLabelsPtr_(nullptr),
75  pointEdgesPtr_(nullptr)
76 {}
77 
78 
79 Foam::faPatch::faPatch
80 (
81  const word& name,
82  const dictionary& dict,
83  const label index,
84  const faBoundaryMesh& bm
85 )
86 :
87  labelList(dict.get<labelList>("edgeLabels")),
88  patchIdentifier(name, dict, index),
89  ngbPolyPatchIndex_(dict.get<label>("ngbPolyPatchIndex")),
90  boundaryMesh_(bm),
91  edgeFacesPtr_(nullptr),
92  pointLabelsPtr_(nullptr),
93  pointEdgesPtr_(nullptr)
94 {}
95 
96 
97 Foam::faPatch::faPatch(const faPatch& p, const faBoundaryMesh& bm)
98 :
99  labelList(p),
100  patchIdentifier(p, p.index()),
101  ngbPolyPatchIndex_(p.ngbPolyPatchIndex_),
102  boundaryMesh_(bm),
103  edgeFacesPtr_(nullptr),
104  pointLabelsPtr_(nullptr),
105  pointEdgesPtr_(nullptr)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
110 
112 {
113  clearOut();
114 }
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
119 Foam::label Foam::faPatch::ngbPolyPatchIndex() const
120 {
121  return ngbPolyPatchIndex_;
122 }
123 
124 
126 {
127  return boundaryMesh_;
128 }
129 
130 
131 Foam::label Foam::faPatch::start() const
132 {
133  return boundaryMesh().mesh().patchStarts()[index()];
134 }
135 
136 
138 {
139  if (!pointLabelsPtr_)
140  {
141  calcPointLabels();
142  }
143 
144  return *pointLabelsPtr_;
145 }
146 
147 
149 {
150  SLList<label> labels;
151 
152  UList<edge> edges = patchSlice(boundaryMesh().mesh().edges());
153 
154  forAll(edges, edgeI)
155  {
156  bool existStart = false;
157  bool existEnd = false;
158 
159  forAllIters(labels, iter)
160  {
161  if (*iter == edges[edgeI].start())
162  {
163  existStart = true;
164  }
165 
166  if (*iter == edges[edgeI].end())
167  {
168  existEnd = true;
169  }
170  }
171 
172  if (!existStart)
173  {
174  labels.append(edges[edgeI].start());
175  }
176 
177  if (!existEnd)
178  {
179  labels.append(edges[edgeI].end());
180  }
181  }
182 
183  pointLabelsPtr_ = new labelList(labels);
184 }
185 
186 
188 {
189  const labelList& points = pointLabels();
190 
191  const edgeList::subList e = patchSlice(boundaryMesh().mesh().edges());
192 
193  // set up storage for pointEdges
194  List<SLList<label>> pointEdgs(points.size());
195 
196  forAll(e, edgeI)
197  {
198  const edge& curPoints = e[edgeI];
199 
200  forAll(curPoints, pointI)
201  {
202  const label localPointIndex = points.find(curPoints[pointI]);
203 
204  pointEdgs[localPointIndex].append(edgeI);
205  }
206  }
207 
208  // sort out the list
209  pointEdgesPtr_ = new labelListList(pointEdgs.size());
210  labelListList& pEdges = *pointEdgesPtr_;
211 
212  forAll(pointEdgs, pointI)
213  {
214  pEdges[pointI].setSize(pointEdgs[pointI].size());
215 
216  label i = 0;
217  for
218  (
219  SLList<label>::iterator curEdgesIter = pointEdgs[pointI].begin();
220  curEdgesIter != pointEdgs[pointI].end();
221  ++curEdgesIter, ++i
222  )
223  {
224  pEdges[pointI][i] = curEdgesIter();
225  }
226  }
227 }
228 
229 
231 {
232  if (!pointEdgesPtr_)
233  {
234  calcPointEdges();
235  }
236 
237  return *pointEdgesPtr_;
238 }
239 
240 
242 {
243  labelList ngbFaces;
244 
245  if (ngbPolyPatchIndex() == -1)
246  {
247  return ngbFaces;
248  }
249 
250  ngbFaces.setSize(faPatch::size());
251 
252  const faMesh& aMesh = boundaryMesh().mesh();
253  const polyMesh& pMesh = aMesh();
255 
256  const labelListList& edgeFaces = pMesh.edgeFaces();
257 
258  labelList faceCells(patch.size(), -1);
259 
260  forAll(faceCells, faceI)
261  {
262  label faceID = aMesh.faceLabels()[faceI];
263 
264  faceCells[faceI] = pMesh.faceOwner()[faceID];
265  }
266 
267  labelList meshEdges
268  (
269  patch.meshEdges
270  (
271  pMesh.edges(),
272  pMesh.cellEdges(),
273  faceCells
274  )
275  );
276 
277  forAll(ngbFaces, edgeI)
278  {
279  ngbFaces[edgeI] = -1;
280 
281  label curEdge = (*this)[edgeI];
282 
283  label curPMeshEdge = meshEdges[curEdge];
284 
285  forAll(edgeFaces[curPMeshEdge], faceI)
286  {
287  label curFace = edgeFaces[curPMeshEdge][faceI];
288 
289  label curPatchID = pMesh.boundaryMesh().whichPatch(curFace);
290 
291  if (curPatchID == ngbPolyPatchIndex())
292  {
293  ngbFaces[edgeI] = curFace;
294  }
295  }
296 
297  if (ngbFaces[edgeI] == -1)
298  {
300  << "Problem with determination of edge ngb faces!"
301  << endl;
302  }
303  }
304 
305  return ngbFaces;
306 }
307 
308 
310 {
311  auto tfN = tmp<vectorField>::New();
312  auto& fN = tfN.ref();
313 
314  if (ngbPolyPatchIndex() == -1)
315  {
316  return tfN;
317  }
318 
319  fN.setSize(faPatch::size());
320 
321  labelList ngbFaces = ngbPolyPatchFaces();
322 
323  const polyMesh& pMesh = boundaryMesh().mesh()();
324 
325  const faceList& faces = pMesh.faces();
326  const pointField& points = pMesh.points();
327 
328  forAll(fN, faceI)
329  {
330  fN[faceI] = faces[ngbFaces[faceI]].unitNormal(points);
331  }
332 
333  return tfN;
334 }
335 
336 
338 {
339  if (ngbPolyPatchIndex() == -1)
340  {
341  return tmp<vectorField>::New();
342  }
343 
344  const labelListList& pntEdges = pointEdges();
345 
346  auto tpN = tmp<vectorField>::New(pntEdges.size(), Zero);
347  auto& pN = tpN.ref();
348 
349  const vectorField faceNormals(ngbPolyPatchFaceNormals());
350 
351  forAll(pN, pointI)
352  {
353  forAll(pntEdges[pointI], edgeI)
354  {
355  pN[pointI] += faceNormals[pntEdges[pointI][edgeI]];
356  }
357  }
358 
359  pN /= mag(pN);
360 
361  return tpN;
362 }
363 
364 
366 {
367  if (!edgeFacesPtr_)
368  {
369  edgeFacesPtr_ = new labelList::subList
370  (
371  patchSlice(boundaryMesh().mesh().edgeOwner())
372  );
373  }
374 
375  return *edgeFacesPtr_;
376 }
377 
378 
380 {
381  return boundaryMesh().mesh().edgeCentres().boundaryField()[index()];
382 }
383 
384 
386 {
387  return boundaryMesh().mesh().Le().boundaryField()[index()];
388 }
389 
390 
392 {
393  return boundaryMesh().mesh().magLe().boundaryField()[index()];
394 }
395 
396 
398 {
399  tmp<vectorField> eN(new vectorField(size()));
400 
401  eN.ref() = edgeLengths()/magEdgeLengths();
402 
403  return eN;
404 }
405 
406 
408 {
409  auto tfc = tmp<vectorField>::New(size());
410  auto& fc = tfc.ref();
411 
412  // get reference to global face centres
413  const vectorField& gfc =
414  boundaryMesh().mesh().areaCentres().internalField();
415 
416  const labelUList& faceLabels = edgeFaces();
417 
418  forAll(faceLabels, edgeI)
419  {
420  fc[edgeI] = gfc[faceLabels[edgeI]];
421  }
422 
423  return tfc;
424 }
425 
426 
428 {
429  return edgeNormals()*(edgeNormals() & (edgeCentres() - edgeFaceCentres()));
430 }
431 
432 
434 {
435  dc = scalar(1)/(edgeNormals() & delta());
436 }
437 
438 
440 {
441  vectorField unitDelta(delta()/mag(delta()));
442  vectorField edgeNormMag(edgeNormals()/mag(edgeNormals()));
443  scalarField dn(edgeNormals() & delta());
444 
445  k = edgeNormMag - (scalar(1)/(unitDelta & edgeNormMag))*unitDelta;
446 }
447 
448 
450 {
451  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
452 }
453 
454 
456 {
457  w = scalar(1);
458 }
459 
460 
462 {
463  return boundaryMesh().mesh().weights().boundaryField()[index()];
464 }
465 
466 
468 {}
469 
470 
471 void Foam::faPatch::resetEdges(const labelList& newEdges)
472 {
473  Info<< "Resetting patch edges" << endl;
474  labelList::operator=(newEdges);
475 
476  clearOut();
477 }
478 
479 
480 void Foam::faPatch::write(Ostream& os) const
481 {
482  os.writeEntry("type", type());
483 
485 
486  const labelList& edgeLabels = *this;
487  edgeLabels.writeEntry("edgeLabels", os);
488  os.writeEntry("ngbPolyPatchIndex", ngbPolyPatchIndex_);
489 }
490 
491 
492 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
493 
495 {
496  p.write(os);
497  os.check(FUNCTION_NAME);
498  return os;
499 }
500 
501 
502 // ************************************************************************* //
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:71
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:206
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::faPatch::ngbPolyPatchIndex
label ngbPolyPatchIndex() const
Return neighbour polyPatch index.
Definition: faPatch.C:120
Foam::faPatch::edgeFaceCentres
tmp< vectorField > edgeFaceCentres() const
Return neighbour face centres.
Definition: faPatch.C:408
Foam::faPatch::pointEdges
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition: faPatch.C:231
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:97
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::faMesh::patch
const indirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMesh.C:915
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:468
Foam::faPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Return patch edge - neighbour face distances.
Definition: faPatch.C:450
Foam::faPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: faPatch.C:456
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
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:149
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:428
Foam::faPatch::weights
const scalarField & weights() const
Return patch weighting factors.
Definition: faPatch.C:462
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:296
Foam::faPatch::ngbPolyPatchFaceNormals
tmp< vectorField > ngbPolyPatchFaceNormals() const
Return normals of neighbour polyPatch faces.
Definition: faPatch.C:310
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:366
Foam::LList
Template class for non-intrusive linked lists.
Definition: LList.H:54
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
Foam::faBoundaryMesh
Finite area boundary mesh.
Definition: faBoundaryMesh.H:66
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::faPatch::boundaryMesh
const faBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: faPatch.C:126
Foam::patchIdentifier
Identifies a patch by name, patch index and physical type.
Definition: patchIdentifier.H:54
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:380
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
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:810
Foam::faPatch::ngbPolyPatchFaces
labelList ngbPolyPatchFaces() const
Return edge neighbour polyPatch faces.
Definition: faPatch.C:242
Foam::List< label >::operator=
void operator=(const UList< label > &a)
Assignment to UList operator. Takes linear time.
Definition: List.C:501
areaFields.H
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::List< label >::subList
SubList< label > subList
Declare type of subList.
Definition: List.H:114
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:121
Foam::faPatch::makeCorrectionVectors
void makeCorrectionVectors(vectorField &) const
Definition: faPatch.C:440
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faPatch::start
label start() const
Patch start in edge list.
Definition: faPatch.C:132
Foam::faPatch::~faPatch
virtual ~faPatch()
Destructor.
Definition: faPatch.C:112
Foam::faPatch::pointLabels
const labelList & pointLabels() const
Return patch point labels.
Definition: faPatch.C:138
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:1094
Foam::faPatch::magEdgeLengths
const scalarField & magEdgeLengths() const
Return edge length magnitudes.
Definition: faPatch.C:392
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:481
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
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
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::faPatch::edgeNormals
tmp< vectorField > edgeNormals() const
Return edge normals.
Definition: faPatch.C:398
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:270
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
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: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::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:338
Foam::faPatch
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:68
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::faPatch::resetEdges
void resetEdges(const labelList &)
Reset edge list.
Definition: faPatch.C:472
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
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::faPatch::edgeLengths
const vectorField & edgeLengths() const
Return edge length vectors.
Definition: faPatch.C:386
Foam::faPatch::calcPointEdges
void calcPointEdges() const
Calculate patch point-edge addressing.
Definition: faPatch.C:188
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:85
Foam::faPatch::size
virtual label size() const
Patch size.
Definition: faPatch.H:249