processorFaPatch.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 "processorFaPatch.H"
30 #include "processorPolyPatch.H" // For newName()
32 #include "transformField.H"
33 #include "faBoundaryMesh.H"
34 #include "faMesh.H"
35 #include "globalMeshData.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
43 }
44 
45 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
46 
48 {
49  neighbPointsPtr_.clear();
50  nonGlobalPatchPointsPtr_.clear();
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const word& name,
59  const labelUList& edgeLabels,
60  const label index,
61  const faBoundaryMesh& bm,
62  const label nbrPolyPatchi,
63  const label myProcNo,
64  const label neighbProcNo,
65  const word& patchType
66 )
67 :
68  coupledFaPatch(name, edgeLabels, index, bm, nbrPolyPatchi, patchType),
69  myProcNo_(myProcNo),
70  neighbProcNo_(neighbProcNo),
71  neighbEdgeCentres_(),
72  neighbEdgeLengths_(),
73  neighbEdgeFaceCentres_(),
74  neighbPointsPtr_(nullptr),
75  nonGlobalPatchPointsPtr_(nullptr)
76 {}
77 
78 
80 (
81  const labelUList& edgeLabels,
82  const label index,
83  const faBoundaryMesh& bm,
84  const label nbrPolyPatchi,
85  const label myProcNo,
86  const label neighbProcNo,
87  const word& patchType
88 )
89 :
91  (
92  processorPolyPatch::newName(myProcNo, neighbProcNo),
93  edgeLabels,
94  index,
95  bm,
96  nbrPolyPatchi,
97  myProcNo,
98  neighbProcNo,
99  patchType
100  )
101 {}
102 
103 
105 (
106  const word& name,
107  const dictionary& dict,
108  const label index,
109  const faBoundaryMesh& bm,
110  const word& patchType
111 )
112 :
113  coupledFaPatch(name, dict, index, bm, patchType),
114  myProcNo_(dict.get<label>("myProcNo")),
115  neighbProcNo_(dict.get<label>("neighbProcNo")),
116  neighbEdgeCentres_(),
117  neighbEdgeLengths_(),
118  neighbEdgeFaceCentres_(),
119  neighbPointsPtr_(nullptr),
120  nonGlobalPatchPointsPtr_(nullptr)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
126 Foam::label Foam::processorFaPatch::comm() const
127 {
128  return boundaryMesh().mesh().comm();
129 }
130 
131 
133 {
134  // If it is not running parallel or there are no global points
135  // create a 1->1 map
136 
137  // Can not use faGlobalMeshData at this point yet
138  // - use polyMesh globalData instead
139 
140  const auto& pMeshGlobalData = boundaryMesh().mesh().mesh().globalData();
141 
142  if (!Pstream::parRun() || !pMeshGlobalData.nGlobalPoints())
143  {
144  // 1 -> 1 mapping
145  nonGlobalPatchPointsPtr_.reset(new labelList(identity(nPoints())));
146  }
147  else
148  {
149  nonGlobalPatchPointsPtr_.reset(new labelList(nPoints()));
150  labelList& ngpp = *nonGlobalPatchPointsPtr_;
151 
152  // Get reference to shared points
153  const labelList& sharedPoints = pMeshGlobalData.sharedPointLabels();
154 
155  const labelList& faMeshPatchPoints = pointLabels();
156 
157  const labelList& meshPoints =
158  boundaryMesh().mesh().patch().meshPoints();
159 
160  label nNonShared = 0;
161 
162  forAll(faMeshPatchPoints, pointi)
163  {
164  const label mpi = meshPoints[faMeshPatchPoints[pointi]];
165  if (!sharedPoints.found(mpi))
166  {
167  ngpp[nNonShared] = pointi;
168  ++nNonShared;
169  }
170  }
171 
172  ngpp.resize(nNonShared);
173  }
174 }
175 
176 
177 void Foam::processorFaPatch::initGeometry(PstreamBuffers& pBufs)
178 {
179  if (Pstream::parRun())
180  {
181  if (neighbProcNo() >= pBufs.nProcs())
182  {
184  << "On patch " << name()
185  << " trying to access out of range neighbour processor "
186  << neighbProcNo() << ". This can happen if" << nl
187  << " trying to run on an incorrect number of processors"
188  << exit(FatalError);
189  }
190 
191  UOPstream toNeighbProc(neighbProcNo(), pBufs);
192 
193  toNeighbProc
195  << edgeLengths()
196  << edgeFaceCentres();
197  }
198 }
199 
200 
201 void Foam::processorFaPatch::calcGeometry(PstreamBuffers& pBufs)
202 {
203  if (Pstream::parRun())
204  {
205  {
206  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
207 
208  fromNeighbProc
209  >> neighbEdgeCentres_
210  >> neighbEdgeLengths_
211  >> neighbEdgeFaceCentres_;
212  }
213 
214  #ifdef FULLDEBUG
215  const scalarField& magEl = magEdgeLengths();
216 
217  forAll(magEl, edgei)
218  {
219  scalar nmagEl = mag(neighbEdgeLengths_[edgei]);
220  scalar avEl = (magEl[edgei] + nmagEl)/2.0;
221 
222  if (mag(magEl[edgei] - nmagEl)/avEl > 1e-6)
223  {
225  << "edge " << edgei
226  << " length does not match neighbour by "
227  << 100*mag(magEl[edgei] - nmagEl)/avEl
228  << "% -- possible edge ordering problem" << nl;
229  }
230  }
231  #endif
232 
233  calcTransformTensors
234  (
235  edgeCentres(),
236  neighbEdgeCentres_,
237  edgeNormals(),
238  neighbEdgeNormals()
239  );
240  }
241 }
242 
243 
245 (
246  PstreamBuffers& pBufs,
247  const pointField& p
248 )
249 {
250  faPatch::movePoints(pBufs, p);
251  initGeometry(pBufs);
252 }
253 
254 
256 (
257  PstreamBuffers& pBufs,
258  const pointField&
259 )
260 {
262 }
263 
264 
266 {
267  // For completeness
269 
270  neighbPointsPtr_.clear();
271 
272  if (Pstream::parRun())
273  {
274  if (neighbProcNo() >= pBufs.nProcs())
275  {
277  << "On patch " << name()
278  << " trying to access out of range neighbour processor "
279  << neighbProcNo() << ". This can happen if" << nl
280  << " trying to run on an incorrect number of processors"
281  << exit(FatalError);
282  }
283 
284  // Express all points as patch edge and index in edge.
285  labelList patchEdge(nPoints());
286  labelList indexInEdge(nPoints());
287 
288  const edgeList::subList patchEdges =
289  patchSlice(boundaryMesh().mesh().edges());
290 
291  const labelListList& ptEdges = pointEdges();
292 
293  for (label patchPointI = 0; patchPointI < nPoints(); ++patchPointI)
294  {
295  label edgeI = ptEdges[patchPointI][0];
296 
297  patchEdge[patchPointI] = edgeI;
298 
299  const edge& e = patchEdges[edgeI];
300 
301  indexInEdge[patchPointI] = e.find(pointLabels()[patchPointI]);
302  }
303 
304  UOPstream toNeighbProc(neighbProcNo(), pBufs);
305 
306  toNeighbProc
307  << patchEdge
308  << indexInEdge;
309  }
310 }
311 
312 
313 void Foam::processorFaPatch::updateMesh(PstreamBuffers& pBufs)
314 {
315  // For completeness
316  faPatch::updateMesh(pBufs);
317 
318  neighbPointsPtr_.clear();
319 
320  if (Pstream::parRun())
321  {
322  labelList nbrPatchEdge;
323  labelList nbrIndexInEdge;
324 
325  {
326  // Note cannot predict exact size since edgeList not (yet) sent as
327  // binary entity but as List of edges.
328 
329  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
330 
331  fromNeighbProc
332  >> nbrPatchEdge
333  >> nbrIndexInEdge;
334  }
335 
336  if (nbrPatchEdge.size() == nPoints())
337  {
338  // Convert neighbour edges and indices into face back into
339  // my edges and points.
340  neighbPointsPtr_.reset(new labelList(nPoints()));
341  labelList& neighbPoints = *neighbPointsPtr_;
342 
343  const edgeList::subList patchEdges =
344  patchSlice(boundaryMesh().mesh().edges());
345 
346  forAll(nbrPatchEdge, nbrPointI)
347  {
348  // Find edge and index in edge on this side.
349  const edge& e = patchEdges[nbrPatchEdge[nbrPointI]];
350 
351  const label index = 1 - nbrIndexInEdge[nbrPointI];
352 
353  const label patchPointI = pointLabels().find(e[index]);
354 
355  neighbPoints[patchPointI] = nbrPointI;
356  }
357  }
358  else
359  {
360  // Differing number of points. Probably patch includes
361  // part of a cyclic.
362  }
363  }
364 }
365 
366 
368 {
369  auto tresult = tmp<vectorField>::New(neighbEdgeLengths_);
370  tresult.ref().normalise();
371  return tresult;
372 }
373 
374 
376 {
377  if (!neighbPointsPtr_)
378  {
379  // Was probably created from cyclic patch and hence the
380  // number of edges or points might differ on both
381  // sides of the processor patch since one side might have
382  // it merged with another bit of geometry
383 
385  << "No extended addressing calculated for patch " << name()
386  << nl
387  << "This can happen if the number of points on both"
388  << " sides of the two coupled patches differ." << nl
389  << "This happens if the processorPatch was constructed from"
390  << " part of a cyclic patch."
392  }
393 
394  return *neighbPointsPtr_;
395 }
396 
397 
399 {
400  if (Pstream::parRun())
401  {
402  // The face normals point in the opposite direction on the other side
403  scalarField neighbEdgeCentresCn
404  (
405  neighbEdgeNormals()
406  & (neighbEdgeCentres() - neighbEdgeFaceCentres())
407  );
408 
409  w = neighbEdgeCentresCn/
410  (
411  (edgeNormals() & coupledFaPatch::delta())
412  + neighbEdgeCentresCn + VSMALL
413  );
414  }
415  else
416  {
417  w = scalar(1);
418  }
419 }
420 
421 
423 {
424  if (Pstream::parRun())
425  {
426  dc = (1.0 - weights())
427  /((edgeNormals() & coupledFaPatch::delta()) + VSMALL);
428  }
429  else
430  {
431  dc = scalar(1)
432  /((edgeNormals() & coupledFaPatch::delta()) + VSMALL);
433  }
434 }
435 
436 
438 {
439  if (Pstream::parRun())
440  {
441  // Do the transformation if necessary
442  if (parallel())
443  {
444  return
446  - (
447  neighbEdgeCentres()
448  - neighbEdgeFaceCentres()
449  );
450  }
451  else
452  {
453  return
455  - transform
456  (
457  forwardT(),
458  (
459  neighbEdgeCentres()
460  - neighbEdgeFaceCentres()
461  )
462  );
463  }
464  }
465  else
466  {
467  return coupledFaPatch::delta();
468  }
469 }
470 
471 
473 {
474  if (!nonGlobalPatchPointsPtr_)
475  {
476  makeNonGlobalPatchPoints();
477  }
478 
479  return *nonGlobalPatchPointsPtr_;
480 }
481 
482 
484 (
485  const labelUList& internalData
486 ) const
487 {
488  return patchInternalField(internalData);
489 }
490 
491 
493 (
494  const labelUList& internalData,
495  const labelUList& edgeFaces
496 ) const
497 {
498  auto tpfld = tmp<labelField>::New();
499  patchInternalField(internalData, edgeFaces, tpfld.ref());
500  return tpfld;
501 }
502 
503 
505 (
506  const Pstream::commsTypes commsType,
507  const labelUList& interfaceData
508 ) const
509 {
510  send(commsType, interfaceData);
511 }
512 
513 
515 (
516  const Pstream::commsTypes commsType,
517  const labelUList&
518 ) const
519 {
520  return receive<label>(commsType, this->size());
521 }
522 
523 
525 (
526  const Pstream::commsTypes commsType,
527  const labelUList& iF
528 ) const
529 {
530  send(commsType, patchInternalField(iF)());
531 }
532 
533 
535 (
536  const Pstream::commsTypes commsType,
537  const labelUList&
538 ) const
539 {
540  return receive<label>(commsType, this->size());
541 }
542 
543 
545 (
546  const Pstream::commsTypes commsType,
547  const labelUList&,
548  const labelUList&
549 ) const
550 {
551  return receive<label>(commsType, this->size());
552 }
553 
554 
556 {
558  os.writeEntry("myProcNo", myProcNo_);
559  os.writeEntry("neighbProcNo", neighbProcNo_);
560 }
561 
562 
563 // ************************************************************************* //
processorFaPatch(const word &name, const labelUList &edgeLabels, const label index, const faBoundaryMesh &bm, const label nbrPolyPatchi, const label myProcNo, const label neighbProcNo, const word &patchType=typeName)
Construct from components with specified name.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList pointLabels(nPoints, -1)
commsTypes
Communications types.
Definition: UPstream.H:72
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
const labelList & nonGlobalPatchPoints() const
Return the set of labels of the processor patch points which are.
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
const labelList & neighbPoints() const
Return neighbour point labels. This is for my local point the.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
SubList< edge > subList
Declare type of subList.
Definition: List.H:144
virtual label comm() const
Return communicator used for communication.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Neighbour processor patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to the interface as a field.
Spatial transformation functions for primitive fields.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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
tmp< vectorField > neighbEdgeNormals() const
Return processor-neighbour patch edge unit normals.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual void write(Ostream &os) const
Write the patch data as a dictionary.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: faPatch.H:174
OBJstream os(runTime.globalPath()/outputName)
void makeNonGlobalPatchPoints() const
Find non-globa patch points.
defineTypeNameAndDebug(combustionModel, 0)
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:72
Buffers for inter-processor communications streams (UOPstream, UIPstream).
#define WarningInFunction
Report a warning using Foam::Warning.
void makeWeights(scalarField &) const
Make patch weighting factors.
coupledFaPatch is an abstract base class for patches that couple regions of the computational domain ...
Finite area boundary mesh.
virtual void initInternalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Initialise neighbour field transfer.
void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
virtual void initTransfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Initialise interface data transfer.
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
virtual tmp< labelField > transfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Transfer and return neighbour field.
const bMesh & mesh() const
Definition: boundaryMesh.H:259
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
Definition: faPatch.C:542
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: faPatch.H:168
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
Processor patch.
virtual ~processorFaPatch()
Destructor.