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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "processorFaPatch.H"
30 #include "IPstream.H"
31 #include "OPstream.H"
32 #include "transformField.H"
33 #include "faBoundaryMesh.H"
34 #include "faMesh.H"
35 #include "globalMeshData.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(processorFaPatch, 0);
42  addToRunTimeSelectionTable(faPatch, processorFaPatch, dictionary);
43 }
44 
45 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
46 
48 {
49  deleteDemandDrivenData(neighbPointsPtr_);
50  deleteDemandDrivenData(nonGlobalPatchPointsPtr_);
51 }
52 
53 
54 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
55 
57 {
58  return boundaryMesh().mesh().comm();
59 }
60 
61 
63 {
64  return Pstream::msgType();
65 }
66 
67 
69 {
70  // If it is not running parallel or there are no global points
71  // create a 1->1 map
72 
73  // Can not use faGlobalMeshData at this point yet
74 
75  if
76  (
78  || !boundaryMesh().mesh()().globalData().nGlobalPoints()
79 // || !boundaryMesh().mesh().globalData().nGlobalPoints()
80  )
81  {
82  nonGlobalPatchPointsPtr_ = new labelList(nPoints());
83  labelList& ngpp = *nonGlobalPatchPointsPtr_;
84  forAll(ngpp, i)
85  {
86  ngpp[i] = i;
87  }
88  }
89  else
90  {
91  // Get reference to shared points
92  const labelList& sharedPoints =
93  boundaryMesh().mesh()().globalData().sharedPointLabels();
94 
95  nonGlobalPatchPointsPtr_ = new labelList(nPoints());
96  labelList& ngpp = *nonGlobalPatchPointsPtr_;
97 
98  const labelList& faMeshPatchPoints = pointLabels();
99 
100  const labelList& meshPoints =
101  boundaryMesh().mesh().patch().meshPoints();
102 
103  label noFiltPoints = 0;
104 
105  forAll(faMeshPatchPoints, pointI)
106  {
107  label curP = meshPoints[faMeshPatchPoints[pointI]];
108 
109  bool found = false;
110 
111  forAll(sharedPoints, sharedI)
112  {
113  if (sharedPoints[sharedI] == curP)
114  {
115  found = true;
116  break;
117  }
118  }
119 
120  if (!found)
121  {
122  ngpp[noFiltPoints] = pointI;
123  ++noFiltPoints;
124  }
125  }
126 
127  ngpp.setSize(noFiltPoints);
128 
129 
130 // // Get reference to shared points
131 // const labelList& sharedPoints =
132 // boundaryMesh().mesh().globalData().sharedPointLabels();
133 
134 // nonGlobalPatchPointsPtr_ = new labelList(nPoints());
135 // labelList& ngpp = *nonGlobalPatchPointsPtr_;
136 
137 // const labelList& patchPoints = pointLabels();
138 
139 // label noFiltPoints = 0;
140 
141 // forAll(patchPoints, pointI)
142 // {
143 // label curP = patchPoints[pointI];
144 
145 // bool found = false;
146 
147 // forAll(sharedPoints, pI)
148 // {
149 // if (sharedPoints[pI] == curP)
150 // {
151 // found = true;
152 // break;
153 // }
154 // }
155 
156 // if (!found)
157 // {
158 // ngpp[noFiltPoints] = pointI;
159 // noFiltPoints++;
160 // }
161 // }
162 
163 // ngpp.setSize(noFiltPoints);
164  }
165 }
166 
167 
169 {
170  if (Pstream::parRun())
171  {
172  OPstream toNeighbProc
173  (
175  neighbProcNo(),
176  3*(sizeof(label) + size()*sizeof(vector))
177  );
178 
179  toNeighbProc
180  << edgeCentres()
181  << edgeLengths()
182  << edgeFaceCentres();
183  }
184 }
185 
186 
188 {
189  if (Pstream::parRun())
190  {
191  {
192  IPstream fromNeighbProc
193  (
195  neighbProcNo(),
196  3*(sizeof(label) + size()*sizeof(vector))
197  );
198  fromNeighbProc
199  >> neighbEdgeCentres_
200  >> neighbEdgeLengths_
201  >> neighbEdgeFaceCentres_;
202  }
203 
204  const scalarField& magEl = magEdgeLengths();
205 
206  forAll(magEl, edgei)
207  {
208  scalar nmagEl = mag(neighbEdgeLengths_[edgei]);
209  scalar avEl = (magEl[edgei] + nmagEl)/2.0;
210 
211  if (mag(magEl[edgei] - nmagEl)/avEl > 1e-6)
212  {
214  << "edge " << edgei
215  << " length does not match neighbour by "
216  << 100*mag(magEl[edgei] - nmagEl)/avEl
217  << "% -- possible edge ordering problem"
218  << exit(FatalError);
219  }
220  }
221 
222  calcTransformTensors
223  (
224  edgeCentres(),
225  neighbEdgeCentres_,
226  edgeNormals(),
227  neighbEdgeLengths_/mag(neighbEdgeLengths_)
228  );
229  }
230 }
231 
232 
234 {
236  initGeometry();
237 }
238 
239 
241 {
242  calcGeometry();
243 }
244 
245 
247 {
248  // For completeness
250 
251  deleteDemandDrivenData(neighbPointsPtr_);
252 
253  if (Pstream::parRun())
254  {
255  // Express all points as patch edge and index in edge.
256  labelList patchEdge(nPoints());
257  labelList indexInEdge(nPoints());
258 
259  const edgeList::subList patchEdges =
260  patchSlice(boundaryMesh().mesh().edges());
261 
262  const labelListList& ptEdges = pointEdges();
263 
264  for (label patchPointI = 0; patchPointI < nPoints(); ++patchPointI)
265  {
266  label edgeI = ptEdges[patchPointI][0];
267 
268  patchEdge[patchPointI] = edgeI;
269 
270  const edge& e = patchEdges[edgeI];
271 
272  indexInEdge[patchPointI] = e.find(pointLabels()[patchPointI]);
273  }
274 
275  OPstream toNeighbProc
276  (
278  neighbProcNo(),
279  2*sizeof(label) + 2*nPoints()*sizeof(label)
280  );
281 
282  toNeighbProc
283  << patchEdge
284  << indexInEdge;
285  }
286 }
287 
288 
290 {
291  // For completeness
293 
294  if (Pstream::parRun())
295  {
296  labelList nbrPatchEdge(nPoints());
297  labelList nbrIndexInEdge(nPoints());
298 
299  {
300  // Note cannot predict exact size since edgeList not (yet) sent as
301  // binary entity but as List of edges.
302  IPstream fromNeighbProc
303  (
305  neighbProcNo()
306  );
307 
308  fromNeighbProc
309  >> nbrPatchEdge
310  >> nbrIndexInEdge;
311  }
312 
313  if (nbrPatchEdge.size() == nPoints())
314  {
315  // Convert neighbour edges and indices into face back into
316  // my edges and points.
317  neighbPointsPtr_ = new labelList(nPoints());
318  labelList& neighbPoints = *neighbPointsPtr_;
319 
320  const edgeList::subList patchEdges =
321  patchSlice(boundaryMesh().mesh().edges());
322 
323  forAll(nbrPatchEdge, nbrPointI)
324  {
325  // Find edge and index in edge on this side.
326  const edge& e = patchEdges[nbrPatchEdge[nbrPointI]];
327 
328  const label index = 1 - nbrIndexInEdge[nbrPointI];
329 
330  const label patchPointI = pointLabels().find(e[index]);
331 
332  neighbPoints[patchPointI] = nbrPointI;
333  }
334  }
335  else
336  {
337  // Differing number of points. Probably patch includes
338  // part of a cyclic.
339  neighbPointsPtr_ = nullptr;
340  }
341  }
342 }
343 
344 
346 {
347  if (!neighbPointsPtr_)
348  {
349  // Was probably created from cyclic patch and hence the
350  // number of edges or points might differ on both
351  // sides of the processor patch since one side might have
352  // it merged with another bit of geometry
353 
355  << "No extended addressing calculated for patch " << name()
356  << nl
357  << "This can happen if the number of points on both"
358  << " sides of the two coupled patches differ." << nl
359  << "This happens if the processorPatch was constructed from"
360  << " part of a cyclic patch."
361  << abort(FatalError);
362  }
363 
364  return *neighbPointsPtr_;
365 }
366 
367 
369 {
370  if (Pstream::parRun())
371  {
372  // The face normals point in the opposite direction on the other side
373  scalarField neighbEdgeCentresCn
374  (
375  (
376  neighbEdgeLengths()
377  /mag(neighbEdgeLengths())
378  )
379  & (
380  neighbEdgeCentres()
381  - neighbEdgeFaceCentres()
382  )
383  );
384 
385  w = neighbEdgeCentresCn/
386  (
387  (edgeNormals() & faPatch::delta())
388  + neighbEdgeCentresCn
389  );
390  }
391  else
392  {
393  w = 1.0;
394  }
395 }
396 
397 
399 {
400  if (Pstream::parRun())
401  {
402  dc = (1.0 - weights())/(edgeNormals() & faPatch::delta());
403  }
404  else
405  {
406  dc = 1.0/(edgeNormals() & faPatch::delta());
407  }
408 }
409 
410 
412 {
413  if (Pstream::parRun())
414  {
415  // To the transformation if necessary
416  if (parallel())
417  {
418  return
420  - (
421  neighbEdgeCentres()
422  - neighbEdgeFaceCentres()
423  );
424  }
425  else
426  {
427  return
429  - transform
430  (
431  forwardT(),
432  (
433  neighbEdgeCentres()
434  - neighbEdgeFaceCentres()
435  )
436  );
437  }
438  }
439  else
440  {
441  return faPatch::delta();
442  }
443 }
444 
445 
447 {
448  if (!nonGlobalPatchPointsPtr_)
449  {
450  makeNonGlobalPatchPoints();
451  }
452 
453  return *nonGlobalPatchPointsPtr_;
454 }
455 
456 
458 (
459  const labelUList& internalData
460 ) const
461 {
462  return patchInternalField(internalData);
463 }
464 
465 
467 (
468  const Pstream::commsTypes commsType,
469  const labelUList& interfaceData
470 ) const
471 {
472  send(commsType, interfaceData);
473 }
474 
475 
477 (
478  const Pstream::commsTypes commsType,
479  const labelUList&
480 ) const
481 {
482  return receive<label>(commsType, this->size());
483 }
484 
485 
487 (
488  const Pstream::commsTypes commsType,
489  const labelUList& iF
490 ) const
491 {
492  send(commsType, patchInternalField(iF)());
493 }
494 
495 
497 (
498  const Pstream::commsTypes commsType,
499  const labelUList&
500 ) const
501 {
502  return receive<label>(commsType, this->size());
503 }
504 
505 
507 {
508  faPatch::write(os);
509  os.writeEntry("myProcNo", myProcNo_);
510  os.writeEntry("neighbProcNo", neighbProcNo_);
511 }
512 
513 
514 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::processorFaPatch::comm
virtual label comm() const
Return communicator used for communication.
Definition: processorFaPatch.C:56
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::UPstream::commsTypes::blocking
Foam::processorFaPatch::movePoints
void movePoints(const pointField &)
Correct patches after moving points.
Definition: processorFaPatch.C:240
Foam::processorFaPatch::updateMesh
virtual void updateMesh()
Update of the patch topology.
Definition: processorFaPatch.C:289
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::processorFaPatch::initGeometry
void initGeometry()
Initialise the calculation of the patch geometry.
Definition: processorFaPatch.C:168
Foam::processorFaPatch::initTransfer
virtual void initTransfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Initialise interface data transfer.
Definition: processorFaPatch.C:467
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::faPatch::movePoints
virtual void movePoints(const pointField &)
Correct patch after moving points.
Definition: faPatch.C:458
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
globalMeshData.H
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
processorFaPatch.H
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
OPstream.H
faMesh.H
Foam::processorFaPatch::transfer
virtual tmp< labelField > transfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Transfer and return neighbour field.
Definition: processorFaPatch.C:477
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::processorFaPatch::initInternalFieldTransfer
virtual void initInternalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Initialise neighbour field transfer.
Definition: processorFaPatch.C:487
Foam::faPatch::delta
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:427
Foam::faPatch::initUpdateMesh
virtual void initUpdateMesh()
Initialise the update of the patch topology.
Definition: faPatch.H:133
transformField.H
Spatial transformation functions for primitive fields.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
IPstream.H
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
Foam::Field< scalar >
Foam::processorFaPatch::neighbPoints
const labelList & neighbPoints() const
Return neighbour point labels. This is for my local point the.
Definition: processorFaPatch.C:345
Foam::processorFaPatch::initUpdateMesh
virtual void initUpdateMesh()
Initialise the update of the patch topology.
Definition: processorFaPatch.C:246
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::FatalError
error FatalError
faBoundaryMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::processorFaPatch::makeNonGlobalPatchPoints
void makeNonGlobalPatchPoints() const
Find non-globa patch points.
Definition: processorFaPatch.C:68
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::faPatch::updateMesh
virtual void updateMesh()
Update of the patch topology.
Definition: faPatch.H:137
Foam::processorFaPatch::write
virtual void write(Ostream &) const
Write the patch data as a dictionary.
Definition: processorFaPatch.C:506
Foam::processorFaPatch::makeDeltaCoeffs
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
Definition: processorFaPatch.C:398
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::processorFaPatch::nonGlobalPatchPoints
const labelList & nonGlobalPatchPoints() const
Return the set of labels of the processor patch points which are.
Definition: processorFaPatch.C:446
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:491
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:66
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::processorFaPatch::~processorFaPatch
virtual ~processorFaPatch()
Destructor.
Definition: processorFaPatch.C:47
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::faPatch::write
virtual void write(Ostream &) const
Write.
Definition: faPatch.C:471
Foam::Vector< scalar >
Foam::processorFaPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to.
Definition: processorFaPatch.C:458
Foam::List< label >
Foam::processorFaPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: processorFaPatch.C:368
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::processorFaPatch::internalFieldTransfer
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
Definition: processorFaPatch.C:497
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
Foam::processorFaPatch::delta
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
Definition: processorFaPatch.C:411
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::processorFaPatch::calcGeometry
void calcGeometry()
Calculate the patch geometry.
Definition: processorFaPatch.C:187
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:418
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::processorFaPatch::initMovePoints
void initMovePoints(const pointField &)
Initialise the patches for moving points.
Definition: processorFaPatch.C:233
pointLabels
labelList pointLabels(nPoints, -1)
Foam::processorFaPatch::tag
virtual int tag() const
Return message tag to use for communication.
Definition: processorFaPatch.C:62