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 #include "demandDrivenData.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(processorFaPatch, 0);
43  addToRunTimeSelectionTable(faPatch, processorFaPatch, dictionary);
44 }
45 
46 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
47 
49 {
50  deleteDemandDrivenData(neighbPointsPtr_);
51  deleteDemandDrivenData(nonGlobalPatchPointsPtr_);
52 }
53 
54 
55 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
56 
57 Foam::label Foam::processorFaPatch::comm() const
58 {
59  return boundaryMesh().mesh().comm();
60 }
61 
62 
64 {
65  return Pstream::msgType();
66 }
67 
68 
70 {
71  // If it is not running parallel or there are no global points
72  // create a 1->1 map
73 
74  // Can not use faGlobalMeshData at this point yet
75 
76  if
77  (
79  || !boundaryMesh().mesh()().globalData().nGlobalPoints()
80 // || !boundaryMesh().mesh().globalData().nGlobalPoints()
81  )
82  {
83  nonGlobalPatchPointsPtr_ = new labelList(nPoints());
84  labelList& ngpp = *nonGlobalPatchPointsPtr_;
85  forAll(ngpp, i)
86  {
87  ngpp[i] = i;
88  }
89  }
90  else
91  {
92  // Get reference to shared points
93  const labelList& sharedPoints =
94  boundaryMesh().mesh()().globalData().sharedPointLabels();
95 
96  nonGlobalPatchPointsPtr_ = new labelList(nPoints());
97  labelList& ngpp = *nonGlobalPatchPointsPtr_;
98 
99  const labelList& faMeshPatchPoints = pointLabels();
100 
101  const labelList& meshPoints =
102  boundaryMesh().mesh().patch().meshPoints();
103 
104  label noFiltPoints = 0;
105 
106  forAll(faMeshPatchPoints, pointI)
107  {
108  label curP = meshPoints[faMeshPatchPoints[pointI]];
109 
110  bool found = false;
111 
112  forAll(sharedPoints, sharedI)
113  {
114  if (sharedPoints[sharedI] == curP)
115  {
116  found = true;
117  break;
118  }
119  }
120 
121  if (!found)
122  {
123  ngpp[noFiltPoints] = pointI;
124  ++noFiltPoints;
125  }
126  }
127 
128  ngpp.setSize(noFiltPoints);
129 
130 
131 // // Get reference to shared points
132 // const labelList& sharedPoints =
133 // boundaryMesh().mesh().globalData().sharedPointLabels();
134 
135 // nonGlobalPatchPointsPtr_ = new labelList(nPoints());
136 // labelList& ngpp = *nonGlobalPatchPointsPtr_;
137 
138 // const labelList& patchPoints = pointLabels();
139 
140 // label noFiltPoints = 0;
141 
142 // forAll(patchPoints, pointI)
143 // {
144 // label curP = patchPoints[pointI];
145 
146 // bool found = false;
147 
148 // forAll(sharedPoints, pI)
149 // {
150 // if (sharedPoints[pI] == curP)
151 // {
152 // found = true;
153 // break;
154 // }
155 // }
156 
157 // if (!found)
158 // {
159 // ngpp[noFiltPoints] = pointI;
160 // noFiltPoints++;
161 // }
162 // }
163 
164 // ngpp.setSize(noFiltPoints);
165  }
166 }
167 
168 
170 {
171  if (Pstream::parRun())
172  {
173  OPstream toNeighbProc
174  (
176  neighbProcNo(),
177  3*(sizeof(label) + size()*sizeof(vector))
178  );
179 
180  toNeighbProc
181  << edgeCentres()
182  << edgeLengths()
183  << edgeFaceCentres();
184  }
185 }
186 
187 
189 {
190  if (Pstream::parRun())
191  {
192  {
193  IPstream fromNeighbProc
194  (
196  neighbProcNo(),
197  3*(sizeof(label) + size()*sizeof(vector))
198  );
199  fromNeighbProc
200  >> neighbEdgeCentres_
201  >> neighbEdgeLengths_
202  >> neighbEdgeFaceCentres_;
203  }
204 
205  const scalarField& magEl = magEdgeLengths();
206 
207  forAll(magEl, edgei)
208  {
209  scalar nmagEl = mag(neighbEdgeLengths_[edgei]);
210  scalar avEl = (magEl[edgei] + nmagEl)/2.0;
211 
212  if (mag(magEl[edgei] - nmagEl)/avEl > 1e-6)
213  {
215  << "edge " << edgei
216  << " length does not match neighbour by "
217  << 100*mag(magEl[edgei] - nmagEl)/avEl
218  << "% -- possible edge ordering problem"
219  << exit(FatalError);
220  }
221  }
222 
223  calcTransformTensors
224  (
225  edgeCentres(),
226  neighbEdgeCentres_,
227  edgeNormals(),
228  neighbEdgeLengths_/mag(neighbEdgeLengths_)
229  );
230  }
231 }
232 
233 
235 {
237  initGeometry();
238 }
239 
240 
242 {
243  calcGeometry();
244 }
245 
246 
248 {
249  // For completeness
251 
252  deleteDemandDrivenData(neighbPointsPtr_);
253 
254  if (Pstream::parRun())
255  {
256  // Express all points as patch edge and index in edge.
257  labelList patchEdge(nPoints());
258  labelList indexInEdge(nPoints());
259 
260  const edgeList::subList patchEdges =
261  patchSlice(boundaryMesh().mesh().edges());
262 
263  const labelListList& ptEdges = pointEdges();
264 
265  for (label patchPointI = 0; patchPointI < nPoints(); ++patchPointI)
266  {
267  label edgeI = ptEdges[patchPointI][0];
268 
269  patchEdge[patchPointI] = edgeI;
270 
271  const edge& e = patchEdges[edgeI];
272 
273  indexInEdge[patchPointI] = e.find(pointLabels()[patchPointI]);
274  }
275 
276  OPstream toNeighbProc
277  (
279  neighbProcNo(),
280  2*sizeof(label) + 2*nPoints()*sizeof(label)
281  );
282 
283  toNeighbProc
284  << patchEdge
285  << indexInEdge;
286  }
287 }
288 
289 
291 {
292  // For completeness
294 
295  if (Pstream::parRun())
296  {
297  labelList nbrPatchEdge(nPoints());
298  labelList nbrIndexInEdge(nPoints());
299 
300  {
301  // Note cannot predict exact size since edgeList not (yet) sent as
302  // binary entity but as List of edges.
303  IPstream fromNeighbProc
304  (
306  neighbProcNo()
307  );
308 
309  fromNeighbProc
310  >> nbrPatchEdge
311  >> nbrIndexInEdge;
312  }
313 
314  if (nbrPatchEdge.size() == nPoints())
315  {
316  // Convert neighbour edges and indices into face back into
317  // my edges and points.
318  neighbPointsPtr_ = new labelList(nPoints());
319  labelList& neighbPoints = *neighbPointsPtr_;
320 
321  const edgeList::subList patchEdges =
322  patchSlice(boundaryMesh().mesh().edges());
323 
324  forAll(nbrPatchEdge, nbrPointI)
325  {
326  // Find edge and index in edge on this side.
327  const edge& e = patchEdges[nbrPatchEdge[nbrPointI]];
328 
329  const label index = 1 - nbrIndexInEdge[nbrPointI];
330 
331  const label patchPointI = pointLabels().find(e[index]);
332 
333  neighbPoints[patchPointI] = nbrPointI;
334  }
335  }
336  else
337  {
338  // Differing number of points. Probably patch includes
339  // part of a cyclic.
340  neighbPointsPtr_ = nullptr;
341  }
342  }
343 }
344 
345 
347 {
348  if (!neighbPointsPtr_)
349  {
350  // Was probably created from cyclic patch and hence the
351  // number of edges or points might differ on both
352  // sides of the processor patch since one side might have
353  // it merged with another bit of geometry
354 
356  << "No extended addressing calculated for patch " << name()
357  << nl
358  << "This can happen if the number of points on both"
359  << " sides of the two coupled patches differ." << nl
360  << "This happens if the processorPatch was constructed from"
361  << " part of a cyclic patch."
362  << abort(FatalError);
363  }
364 
365  return *neighbPointsPtr_;
366 }
367 
368 
370 {
371  if (Pstream::parRun())
372  {
373  // The face normals point in the opposite direction on the other side
374  scalarField neighbEdgeCentresCn
375  (
376  (
377  neighbEdgeLengths()
378  /mag(neighbEdgeLengths())
379  )
380  & (
381  neighbEdgeCentres()
382  - neighbEdgeFaceCentres()
383  )
384  );
385 
386  w = neighbEdgeCentresCn/
387  (
388  (edgeNormals() & faPatch::delta())
389  + neighbEdgeCentresCn
390  );
391  }
392  else
393  {
394  w = 1.0;
395  }
396 }
397 
398 
400 {
401  if (Pstream::parRun())
402  {
403  dc = (1.0 - weights())/(edgeNormals() & faPatch::delta());
404  }
405  else
406  {
407  dc = 1.0/(edgeNormals() & faPatch::delta());
408  }
409 }
410 
411 
413 {
414  if (Pstream::parRun())
415  {
416  // To the transformation if necessary
417  if (parallel())
418  {
419  return
421  - (
422  neighbEdgeCentres()
423  - neighbEdgeFaceCentres()
424  );
425  }
426  else
427  {
428  return
430  - transform
431  (
432  forwardT(),
433  (
434  neighbEdgeCentres()
435  - neighbEdgeFaceCentres()
436  )
437  );
438  }
439  }
440  else
441  {
442  return faPatch::delta();
443  }
444 }
445 
446 
448 {
449  if (!nonGlobalPatchPointsPtr_)
450  {
451  makeNonGlobalPatchPoints();
452  }
453 
454  return *nonGlobalPatchPointsPtr_;
455 }
456 
457 
459 (
460  const labelUList& internalData
461 ) const
462 {
463  return patchInternalField(internalData);
464 }
465 
466 
468 (
469  const Pstream::commsTypes commsType,
470  const labelUList& interfaceData
471 ) const
472 {
473  send(commsType, interfaceData);
474 }
475 
476 
478 (
479  const Pstream::commsTypes commsType,
480  const labelUList&
481 ) const
482 {
483  return receive<label>(commsType, this->size());
484 }
485 
486 
488 (
489  const Pstream::commsTypes commsType,
490  const labelUList& iF
491 ) const
492 {
493  send(commsType, patchInternalField(iF)());
494 }
495 
496 
498 (
499  const Pstream::commsTypes commsType,
500  const labelUList&
501 ) const
502 {
503  return receive<label>(commsType, this->size());
504 }
505 
506 
508 {
509  faPatch::write(os);
510  os.writeEntry("myProcNo", myProcNo_);
511  os.writeEntry("neighbProcNo", neighbProcNo_);
512 }
513 
514 
515 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::processorFaPatch::comm
virtual label comm() const
Return communicator used for communication.
Definition: processorFaPatch.C:57
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::UPstream::commsTypes::blocking
Foam::processorFaPatch::movePoints
void movePoints(const pointField &)
Correct patches after moving points.
Definition: processorFaPatch.C:241
Foam::processorFaPatch::updateMesh
virtual void updateMesh()
Update of the patch topology.
Definition: processorFaPatch.C:290
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::processorFaPatch::initGeometry
void initGeometry()
Initialise the calculation of the patch geometry.
Definition: processorFaPatch.C:169
Foam::processorFaPatch::initTransfer
virtual void initTransfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Initialise interface data transfer.
Definition: processorFaPatch.C:468
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::faPatch::movePoints
virtual void movePoints(const pointField &)
Correct patch after moving points.
Definition: faPatch.C:468
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
globalMeshData.H
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
processorFaPatch.H
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
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:478
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:488
Foam::faPatch::delta
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:428
Foam::faPatch::initUpdateMesh
virtual void initUpdateMesh()
Initialise the update of the patch topology.
Definition: faPatch.H:136
transformField.H
Spatial transformation functions for primitive fields.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
IPstream.H
Foam::Field< scalar >
Foam::processorFaPatch::neighbPoints
const labelList & neighbPoints() const
Return neighbour point labels. This is for my local point the.
Definition: processorFaPatch.C:346
Foam::processorFaPatch::initUpdateMesh
virtual void initUpdateMesh()
Initialise the update of the patch topology.
Definition: processorFaPatch.C:247
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:144
Foam::processorFaPatch::makeNonGlobalPatchPoints
void makeNonGlobalPatchPoints() const
Find non-globa patch points.
Definition: processorFaPatch.C:69
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:140
Foam::processorFaPatch::write
virtual void write(Ostream &) const
Write the patch data as a dictionary.
Definition: processorFaPatch.C:507
Foam::processorFaPatch::makeDeltaCoeffs
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
Definition: processorFaPatch.C:399
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:447
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:541
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::processorFaPatch::~processorFaPatch
virtual ~processorFaPatch()
Destructor.
Definition: processorFaPatch.C:48
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::faPatch::write
virtual void write(Ostream &) const
Write.
Definition: faPatch.C:481
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:459
Foam::List< label >
Foam::processorFaPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: processorFaPatch.C:369
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:498
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:232
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
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:412
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:188
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:310
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:234
pointLabels
labelList pointLabels(nPoints, -1)
Foam::processorFaPatch::tag
virtual int tag() const
Return message tag to use for communication.
Definition: processorFaPatch.C:63