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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
39namespace 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 * * * * * * * * * * * * * * //
125
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
139 if
140 (
142 || !boundaryMesh().mesh()().globalData().nGlobalPoints()
143 )
144 {
145 // 1 -> 1 mapping
146 nonGlobalPatchPointsPtr_.reset(new labelList(identity(nPoints())));
147 }
148 else
149 {
150 nonGlobalPatchPointsPtr_.reset(new labelList(nPoints()));
151 labelList& ngpp = *nonGlobalPatchPointsPtr_;
152
153 // Get reference to shared points
154 const labelList& sharedPoints =
155 boundaryMesh().mesh()().globalData().sharedPointLabels();
156
157 const labelList& faMeshPatchPoints = pointLabels();
158
159 const labelList& meshPoints =
160 boundaryMesh().mesh().patch().meshPoints();
161
162 label nNonShared = 0;
163
164 forAll(faMeshPatchPoints, pointi)
165 {
166 const label mpi = meshPoints[faMeshPatchPoints[pointi]];
167 if (!sharedPoints.found(mpi))
168 {
169 ngpp[nNonShared] = pointi;
170 ++nNonShared;
171 }
172 }
173
174 ngpp.resize(nNonShared);
175 }
176}
177
178
180{
181 if (Pstream::parRun())
182 {
183 if (neighbProcNo() >= Pstream::nProcs(pBufs.comm()))
184 {
186 << "On patch " << name()
187 << " trying to access out of range neighbour processor "
188 << neighbProcNo() << ". This can happen if" << nl
189 << " trying to run on an incorrect number of processors"
190 << exit(FatalError);
191 }
192
193 UOPstream toNeighbProc(neighbProcNo(), pBufs);
194
195 toNeighbProc
196 << edgeCentres()
197 << edgeLengths()
198 << edgeFaceCentres();
199 }
200}
201
202
204{
205 if (Pstream::parRun())
206 {
207 {
208 UIPstream fromNeighbProc(neighbProcNo(), pBufs);
209
210 fromNeighbProc
211 >> neighbEdgeCentres_
212 >> neighbEdgeLengths_
213 >> neighbEdgeFaceCentres_;
214 }
215
216 const scalarField& magEl = magEdgeLengths();
217
218 forAll(magEl, edgei)
219 {
220 scalar nmagEl = mag(neighbEdgeLengths_[edgei]);
221 scalar avEl = (magEl[edgei] + nmagEl)/2.0;
222
223 if (mag(magEl[edgei] - nmagEl)/avEl > 1e-6)
224 {
226 << "edge " << edgei
227 << " length does not match neighbour by "
228 << 100*mag(magEl[edgei] - nmagEl)/avEl
229 << "% -- possible edge ordering problem"
230 << exit(FatalError);
231 }
232 }
233
234 calcTransformTensors
235 (
236 edgeCentres(),
237 neighbEdgeCentres_,
238 edgeNormals(),
239 neighbEdgeLengths_/mag(neighbEdgeLengths_)
240 );
241 }
242}
243
244
246(
247 PstreamBuffers& pBufs,
248 const pointField& p
249)
250{
251 faPatch::movePoints(pBufs, p);
252 initGeometry(pBufs);
253}
254
255
257(
258 PstreamBuffers& pBufs,
259 const pointField&
260)
261{
263}
264
265
267{
268 // For completeness
270
271 neighbPointsPtr_.clear();
272
273 if (Pstream::parRun())
274 {
275 if (neighbProcNo() >= Pstream::nProcs(pBufs.comm()))
276 {
278 << "On patch " << name()
279 << " trying to access out of range neighbour processor "
280 << neighbProcNo() << ". This can happen if" << nl
281 << " trying to run on an incorrect number of processors"
282 << exit(FatalError);
283 }
284
285 // Express all points as patch edge and index in edge.
286 labelList patchEdge(nPoints());
287 labelList indexInEdge(nPoints());
288
289 const edgeList::subList patchEdges =
290 patchSlice(boundaryMesh().mesh().edges());
291
292 const labelListList& ptEdges = pointEdges();
293
294 for (label patchPointI = 0; patchPointI < nPoints(); ++patchPointI)
295 {
296 label edgeI = ptEdges[patchPointI][0];
297
298 patchEdge[patchPointI] = edgeI;
299
300 const edge& e = patchEdges[edgeI];
301
302 indexInEdge[patchPointI] = e.find(pointLabels()[patchPointI]);
303 }
304
305 UOPstream toNeighbProc(neighbProcNo(), pBufs);
306
307 toNeighbProc
308 << patchEdge
309 << indexInEdge;
310 }
311}
312
313
315{
316 // For completeness
317 faPatch::updateMesh(pBufs);
318
319 neighbPointsPtr_.clear();
320
321 if (Pstream::parRun())
322 {
323 labelList nbrPatchEdge(nPoints());
324 labelList nbrIndexInEdge(nPoints());
325
326 {
327 // Note cannot predict exact size since edgeList not (yet) sent as
328 // binary entity but as List of edges.
329
330 UIPstream fromNeighbProc(neighbProcNo(), pBufs);
331
332 fromNeighbProc
333 >> nbrPatchEdge
334 >> nbrIndexInEdge;
335 }
336
337 if (nbrPatchEdge.size() == nPoints())
338 {
339 // Convert neighbour edges and indices into face back into
340 // my edges and points.
341 neighbPointsPtr_.reset(new labelList(nPoints()));
342 labelList& neighbPoints = *neighbPointsPtr_;
343
344 const edgeList::subList patchEdges =
345 patchSlice(boundaryMesh().mesh().edges());
346
347 forAll(nbrPatchEdge, nbrPointI)
348 {
349 // Find edge and index in edge on this side.
350 const edge& e = patchEdges[nbrPatchEdge[nbrPointI]];
351
352 const label index = 1 - nbrIndexInEdge[nbrPointI];
353
354 const label patchPointI = pointLabels().find(e[index]);
355
356 neighbPoints[patchPointI] = nbrPointI;
357 }
358 }
359 else
360 {
361 // Differing number of points. Probably patch includes
362 // part of a cyclic.
363 }
364 }
365}
366
367
369{
370 if (!neighbPointsPtr_)
371 {
372 // Was probably created from cyclic patch and hence the
373 // number of edges or points might differ on both
374 // sides of the processor patch since one side might have
375 // it merged with another bit of geometry
376
378 << "No extended addressing calculated for patch " << name()
379 << nl
380 << "This can happen if the number of points on both"
381 << " sides of the two coupled patches differ." << nl
382 << "This happens if the processorPatch was constructed from"
383 << " part of a cyclic patch."
384 << abort(FatalError);
385 }
386
387 return *neighbPointsPtr_;
388}
389
390
392{
393 if (Pstream::parRun())
394 {
395 // The face normals point in the opposite direction on the other side
396 scalarField neighbEdgeCentresCn
397 (
398 (
399 neighbEdgeLengths()
400 /mag(neighbEdgeLengths())
401 )
402 & (
403 neighbEdgeCentres()
404 - neighbEdgeFaceCentres()
405 )
406 );
407
408 w = neighbEdgeCentresCn/
409 (
410 (edgeNormals() & faPatch::delta())
411 + neighbEdgeCentresCn
412 );
413 }
414 else
415 {
416 w = 1.0;
417 }
418}
419
420
422{
423 if (Pstream::parRun())
424 {
425 dc = (1.0 - weights())/(edgeNormals() & faPatch::delta());
426 }
427 else
428 {
429 dc = 1.0/(edgeNormals() & faPatch::delta());
430 }
431}
432
433
435{
436 if (Pstream::parRun())
437 {
438 // To the transformation if necessary
439 if (parallel())
440 {
441 return
443 - (
444 neighbEdgeCentres()
445 - neighbEdgeFaceCentres()
446 );
447 }
448 else
449 {
450 return
452 - transform
453 (
454 forwardT(),
455 (
456 neighbEdgeCentres()
457 - neighbEdgeFaceCentres()
458 )
459 );
460 }
461 }
462 else
463 {
464 return faPatch::delta();
465 }
466}
467
468
470{
471 if (!nonGlobalPatchPointsPtr_)
472 {
473 makeNonGlobalPatchPoints();
474 }
475
476 return *nonGlobalPatchPointsPtr_;
477}
478
479
481(
482 const labelUList& internalData
483) const
484{
485 return patchInternalField(internalData);
486}
487
488
490(
491 const labelUList& internalData,
492 const labelUList& edgeFaces
493) const
494{
495 return patchInternalField(internalData, edgeFaces);
496}
497
498
500(
501 const Pstream::commsTypes commsType,
502 const labelUList& interfaceData
503) const
504{
505 send(commsType, interfaceData);
506}
507
508
510(
511 const Pstream::commsTypes commsType,
512 const labelUList&
513) const
514{
515 return receive<label>(commsType, this->size());
516}
517
518
520(
521 const Pstream::commsTypes commsType,
522 const labelUList& iF
523) const
524{
525 send(commsType, patchInternalField(iF)());
526}
527
528
530(
531 const Pstream::commsTypes commsType,
532 const labelUList&
533) const
534{
535 return receive<label>(commsType, this->size());
536}
537
538
540(
541 const Pstream::commsTypes commsType,
542 const labelUList&,
543 const labelUList&
544) const
545{
546 return receive<label>(commsType, this->size());
547}
548
549
551{
553 os.writeEntry("myProcNo", myProcNo_);
554 os.writeEntry("neighbProcNo", neighbProcNo_);
555}
556
557
558// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
label comm() const noexcept
Communicator.
A List obtained as a section of another List.
Definition: SubList.H:70
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:265
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Types of communications.
Definition: UPstream.H:67
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
const bMesh & mesh() const
Definition: boundaryMesh.H:206
coupledFaPatch is an abstract base class for patches that couple regions of the computational domain ...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Finite area boundary mesh.
void calcGeometry()
Calculate the geometry for the patches.
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:78
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: faPatch.C:494
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: faPatch.H:138
virtual bool write()
Write the output fields.
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:188
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
virtual void initTransfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Initialise interface data transfer.
virtual ~processorFaPatch()
Destructor.
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
void makeWeights(scalarField &) const
Make patch weighting factors.
const labelList & neighbPoints() const
Return neighbour point labels. This is for my local point the.
void makeNonGlobalPatchPoints() const
Find non-globa patch points.
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual label comm() const
Return communicator used for communication.
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
virtual void initInternalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Initialise neighbour field transfer.
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
const labelList & nonGlobalPatchPoints() const
Return the set of labels of the processor patch points which are.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Neighbour processor patch.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
label nPoints
Namespace for OpenFOAM.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
List< label > labelList
A List of labels.
Definition: List.H:66
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList pointLabels(nPoints, -1)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Spatial transformation functions for primitive fields.