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 "faMesh.H"
34 #include "globalMeshData.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
42 }
43 
44 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45 
47 {
48  neighbPointsPtr_.clear();
49  nonGlobalPatchPointsPtr_.clear();
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const word& name,
58  const labelUList& edgeLabels,
59  const label index,
60  const faBoundaryMesh& bm,
61  const label nbrPolyPatchi,
62  const label myProcNo,
63  const label neighbProcNo,
64  const word& patchType
65 )
66 :
67  coupledFaPatch(name, edgeLabels, index, bm, nbrPolyPatchi, patchType),
68  myProcNo_(myProcNo),
69  neighbProcNo_(neighbProcNo),
70  neighbEdgeCentres_(),
71  neighbEdgeLengths_(),
72  neighbEdgeFaceCentres_(),
73  neighbPointsPtr_(nullptr),
74  nonGlobalPatchPointsPtr_(nullptr)
75 {}
76 
77 
79 (
80  const labelUList& edgeLabels,
81  const label index,
82  const faBoundaryMesh& bm,
83  const label nbrPolyPatchi,
84  const label myProcNo,
85  const label neighbProcNo,
86  const word& patchType
87 )
88 :
90  (
91  processorPolyPatch::newName(myProcNo, neighbProcNo),
92  edgeLabels,
93  index,
94  bm,
95  nbrPolyPatchi,
96  myProcNo,
97  neighbProcNo,
98  patchType
99  )
100 {}
101 
102 
104 (
105  const word& name,
106  const dictionary& dict,
107  const label index,
108  const faBoundaryMesh& bm,
109  const word& patchType
110 )
111 :
112  coupledFaPatch(name, dict, index, bm, patchType),
113  myProcNo_(dict.get<label>("myProcNo")),
114  neighbProcNo_(dict.get<label>("neighbProcNo")),
115  neighbEdgeCentres_(),
116  neighbEdgeLengths_(),
117  neighbEdgeFaceCentres_(),
118  neighbPointsPtr_(nullptr),
119  nonGlobalPatchPointsPtr_(nullptr)
120 {}
121 
122 
123 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125 Foam::label Foam::processorFaPatch::comm() const
126 {
127  return boundaryMesh().mesh().comm();
128 }
129 
130 
132 {
133  // If it is not running parallel or there are no global points
134  // create a 1->1 map
135 
136  // Can not use faGlobalMeshData at this point yet
137  // - use polyMesh globalData instead
138 
139  const auto& pMeshGlobalData = boundaryMesh().mesh().mesh().globalData();
140 
141  if (!Pstream::parRun() || !pMeshGlobalData.nGlobalPoints())
142  {
143  // 1 -> 1 mapping
144  nonGlobalPatchPointsPtr_.reset(new labelList(identity(nPoints())));
145  }
146  else
147  {
148  nonGlobalPatchPointsPtr_.reset(new labelList(nPoints()));
149  labelList& ngpp = *nonGlobalPatchPointsPtr_;
150 
151  // Get reference to shared points
152  const labelList& sharedPoints = pMeshGlobalData.sharedPointLabels();
153 
154  const labelList& faMeshPatchPoints = pointLabels();
155 
156  const labelList& meshPoints =
157  boundaryMesh().mesh().patch().meshPoints();
158 
159  label nNonShared = 0;
160 
161  forAll(faMeshPatchPoints, pointi)
162  {
163  const label mpi = meshPoints[faMeshPatchPoints[pointi]];
164  if (!sharedPoints.found(mpi))
165  {
166  ngpp[nNonShared] = pointi;
167  ++nNonShared;
168  }
169  }
170 
171  ngpp.resize(nNonShared);
172  }
173 }
174 
175 
176 void Foam::processorFaPatch::initGeometry(PstreamBuffers& pBufs)
177 {
178  if (Pstream::parRun())
179  {
180  if (neighbProcNo() >= pBufs.nProcs())
181  {
183  << "On patch " << name()
184  << " trying to access out of range neighbour processor "
185  << neighbProcNo() << ". This can happen if" << nl
186  << " trying to run on an incorrect number of processors"
187  << exit(FatalError);
188  }
189 
190  UOPstream toNeighbProc(neighbProcNo(), pBufs);
191 
192  toNeighbProc
194  << edgeLengths()
195  << edgeFaceCentres();
196  }
197 }
198 
199 
200 void Foam::processorFaPatch::calcGeometry(PstreamBuffers& pBufs)
201 {
202  if (Pstream::parRun())
203  {
204  {
205  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
206 
207  fromNeighbProc
208  >> neighbEdgeCentres_
209  >> neighbEdgeLengths_
210  >> neighbEdgeFaceCentres_;
211  }
212 
213  #ifdef FULLDEBUG
214  const scalarField& magEl = magEdgeLengths();
215 
216  forAll(magEl, edgei)
217  {
218  scalar nmagEl = mag(neighbEdgeLengths_[edgei]);
219  scalar avEl = (magEl[edgei] + nmagEl)/2.0;
220 
221  if (mag(magEl[edgei] - nmagEl)/avEl > 1e-6)
222  {
224  << "edge " << edgei
225  << " length does not match neighbour by "
226  << 100*mag(magEl[edgei] - nmagEl)/avEl
227  << "% -- possible edge ordering problem" << nl;
228  }
229  }
230  #endif
231 
232  calcTransformTensors
233  (
234  edgeCentres(),
235  neighbEdgeCentres_,
236  edgeNormals(),
237  neighbEdgeNormals()
238  );
239  }
240 }
241 
242 
244 (
245  PstreamBuffers& pBufs,
246  const pointField& p
247 )
248 {
249  faPatch::movePoints(pBufs, p);
250  initGeometry(pBufs);
251 }
252 
253 
255 (
256  PstreamBuffers& pBufs,
257  const pointField&
258 )
259 {
261 }
262 
263 
265 {
266  // For completeness
268 
269  neighbPointsPtr_.clear();
270 
271  if (Pstream::parRun())
272  {
273  if (neighbProcNo() >= pBufs.nProcs())
274  {
276  << "On patch " << name()
277  << " trying to access out of range neighbour processor "
278  << neighbProcNo() << ". This can happen if" << nl
279  << " trying to run on an incorrect number of processors"
280  << exit(FatalError);
281  }
282 
283  // Express all points as patch edge and index in edge.
284  labelList patchEdge(nPoints());
285  labelList indexInEdge(nPoints());
286 
287  const edgeList::subList patchEdges =
288  patchSlice(boundaryMesh().mesh().edges());
289 
290  const labelListList& ptEdges = pointEdges();
291 
292  for (label patchPointI = 0; patchPointI < nPoints(); ++patchPointI)
293  {
294  label edgeI = ptEdges[patchPointI][0];
295 
296  patchEdge[patchPointI] = edgeI;
297 
298  const edge& e = patchEdges[edgeI];
299 
300  indexInEdge[patchPointI] = e.find(pointLabels()[patchPointI]);
301  }
302 
303  UOPstream toNeighbProc(neighbProcNo(), pBufs);
304 
305  toNeighbProc
306  << patchEdge
307  << indexInEdge;
308  }
309 }
310 
311 
312 void Foam::processorFaPatch::updateMesh(PstreamBuffers& pBufs)
313 {
314  // For completeness
315  faPatch::updateMesh(pBufs);
316 
317  neighbPointsPtr_.clear();
318 
319  if (Pstream::parRun())
320  {
321  labelList nbrPatchEdge;
322  labelList nbrIndexInEdge;
323 
324  {
325  // Note cannot predict exact size since edgeList not (yet) sent as
326  // binary entity but as List of edges.
327 
328  UIPstream fromNeighbProc(neighbProcNo(), pBufs);
329 
330  fromNeighbProc
331  >> nbrPatchEdge
332  >> nbrIndexInEdge;
333  }
334 
335  if (nbrPatchEdge.size() == nPoints())
336  {
337  // Convert neighbour edges and indices into face back into
338  // my edges and points.
339  neighbPointsPtr_.reset(new labelList(nPoints()));
340  labelList& neighbPoints = *neighbPointsPtr_;
341 
342  const edgeList::subList patchEdges =
343  patchSlice(boundaryMesh().mesh().edges());
344 
345  forAll(nbrPatchEdge, nbrPointI)
346  {
347  // Find edge and index in edge on this side.
348  const edge& e = patchEdges[nbrPatchEdge[nbrPointI]];
349 
350  const label index = 1 - nbrIndexInEdge[nbrPointI];
351 
352  const label patchPointI = pointLabels().find(e[index]);
353 
354  neighbPoints[patchPointI] = nbrPointI;
355  }
356  }
357  else
358  {
359  // Differing number of points. Probably patch includes
360  // part of a cyclic.
361  }
362  }
363 }
364 
365 
367 {
368  auto tresult = tmp<vectorField>::New(neighbEdgeLengths_);
369  tresult.ref().normalise();
370  return tresult;
371 }
372 
373 
375 {
376  if (!neighbPointsPtr_)
377  {
378  // Was probably created from cyclic patch and hence the
379  // number of edges or points might differ on both
380  // sides of the processor patch since one side might have
381  // it merged with another bit of geometry
382 
384  << "No extended addressing calculated for patch " << name()
385  << nl
386  << "This can happen if the number of points on both"
387  << " sides of the two coupled patches differ." << nl
388  << "This happens if the processorPatch was constructed from"
389  << " part of a cyclic patch."
391  }
392 
393  return *neighbPointsPtr_;
394 }
395 
396 
398 {
399  if (Pstream::parRun())
400  {
401  // The face normals point in the opposite direction on the other side
402  scalarField neighbEdgeCentresCn
403  (
404  neighbEdgeNormals()
405  & (neighbEdgeCentres() - neighbEdgeFaceCentres())
406  );
407 
408  w = neighbEdgeCentresCn/
409  (
410  (edgeNormals() & coupledFaPatch::delta())
411  + neighbEdgeCentresCn + VSMALL
412  );
413  }
414  else
415  {
416  w = scalar(1);
417  }
418 }
419 
420 
422 {
423  if (Pstream::parRun())
424  {
425  dc = (1.0 - weights())
426  /((edgeNormals() & coupledFaPatch::delta()) + VSMALL);
427  }
428  else
429  {
430  dc = scalar(1)
431  /((edgeNormals() & coupledFaPatch::delta()) + VSMALL);
432  }
433 }
434 
435 
437 {
438  if (Pstream::parRun())
439  {
440  // Do the transformation if necessary
441  if (parallel())
442  {
443  return
445  - (
446  neighbEdgeCentres()
447  - neighbEdgeFaceCentres()
448  );
449  }
450  else
451  {
452  return
454  - transform
455  (
456  forwardT(),
457  (
458  neighbEdgeCentres()
459  - neighbEdgeFaceCentres()
460  )
461  );
462  }
463  }
464  else
465  {
466  return coupledFaPatch::delta();
467  }
468 }
469 
470 
472 {
473  if (!nonGlobalPatchPointsPtr_)
474  {
475  makeNonGlobalPatchPoints();
476  }
477 
478  return *nonGlobalPatchPointsPtr_;
479 }
480 
481 
483 (
484  const labelUList& internalData
485 ) const
486 {
487  return patchInternalField(internalData);
488 }
489 
490 
492 (
493  const labelUList& internalData,
494  const labelUList& edgeFaces
495 ) const
496 {
497  auto tpfld = tmp<labelField>::New();
498  patchInternalField(internalData, edgeFaces, tpfld.ref());
499  return tpfld;
500 }
501 
502 
504 (
505  const Pstream::commsTypes commsType,
506  const labelUList& interfaceData
507 ) const
508 {
509  send(commsType, interfaceData);
510 }
511 
512 
514 (
515  const Pstream::commsTypes commsType,
516  const labelUList&
517 ) const
518 {
519  return receive<label>(commsType, this->size());
520 }
521 
522 
524 (
525  const Pstream::commsTypes commsType,
526  const labelUList& iF
527 ) const
528 {
529  send(commsType, patchInternalField(iF)());
530 }
531 
532 
534 (
535  const Pstream::commsTypes commsType,
536  const labelUList&
537 ) const
538 {
539  return receive<label>(commsType, this->size());
540 }
541 
542 
544 (
545  const Pstream::commsTypes commsType,
546  const labelUList&,
547  const labelUList&
548 ) const
549 {
550  return receive<label>(commsType, this->size());
551 }
552 
553 
555 {
557  os.writeEntry("myProcNo", myProcNo_);
558  os.writeEntry("neighbProcNo", neighbProcNo_);
559 }
560 
561 
562 // ************************************************************************* //
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:77
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:608
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:1061
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:546
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:528
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.