faMeshTopology.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) 2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "faMesh.H"
29#include "faMeshBoundaryHalo.H"
30#include "globalMeshData.H"
32#include "edgeHashes.H"
33#include "foamVtkLineWriter.H"
35
36// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// Print out edges as point pairs
42template<class PatchType>
43static void printPatchEdges
44(
45 Ostream& os,
46 const PatchType& p,
47 const labelList& edgeIds,
48 label maxOutput = 10
49)
50{
51 label nOutput = 0;
52
53 for (const label patchEdgei : edgeIds)
54 {
55 const edge e(p.meshEdge(patchEdgei));
56
57 os << " "
58 << p.points()[e.first()] << ' '
59 << p.points()[e.second()] << nl;
60
61 ++nOutput;
62 if (maxOutput > 0 && nOutput >= maxOutput)
63 {
64 os << " ... suppressing further output" << nl;
65 break;
66 }
67 }
68}
69
70} // End namespace Foam
71
72
73// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
74
76Foam::faMesh::getBoundaryEdgeConnections() const
77{
78 const polyBoundaryMesh& pbm = mesh().boundaryMesh();
79
80 const label nNonProcessor = pbm.nNonProcessor();
81
82 const label nInternalEdges = patch().nInternalEdges();
83 const label nBoundaryEdges = patch().nBoundaryEdges();
84
85 // The output result:
86 List<Pair<patchTuple>> bndEdgeConnections(nBoundaryEdges);
87
88 // Map edges (mesh numbering) back to a boundary index
89 EdgeMap<label> edgeToBoundaryIndex(2*nBoundaryEdges);
90
91 label nBadEdges(0);
92 labelHashSet badEdges(2*nBoundaryEdges);
93
94 {
95 // Local collection structure for accounting of patch pairs.
96 // Based on 'edge' which has some hash-like insertion properties
97 // that are useful.
98 struct patchPairingType : public Foam::edge
99 {
100 label patchEdgei_ = -1;
101 label meshFacei_ = -1;
102
103 void clear()
104 {
106 patchEdgei_ = -1;
107 meshFacei_ = -1;
108 }
109 };
110
111 List<patchPairingType> patchPairings(nBoundaryEdges);
112
114 << "Determining required boundary edge connections, "
115 << "resolving locally attached boundary edges." << endl;
116
117
118 // Pass 1:
119 // - setup lookup (edge -> bnd index)
120 // - add owner patch for each boundary edge
121 for (label bndEdgei = 0; bndEdgei < nBoundaryEdges; ++bndEdgei)
122 {
123 const label patchEdgei = (bndEdgei + nInternalEdges);
124
125 edgeToBoundaryIndex.insert
126 (
127 patch().meshEdge(patchEdgei),
128 bndEdgei
129 );
130
131 // The attached patch face. Should only be one!
132 const labelList& edgeFaces = patch().edgeFaces()[patchEdgei];
133
134 if (edgeFaces.size() != 1)
135 {
136 badEdges.insert(patchEdgei);
137 continue;
138 }
139
140 const label patchFacei = edgeFaces[0];
141 const label meshFacei = faceLabels_[patchFacei];
142 const label bndFacei = (meshFacei - mesh().nInternalFaces());
143
145 const label patchId = pbm.patchID()[bndFacei];
146
147 // Primary bookkeeping
148 {
149 auto& tuple = bndEdgeConnections[bndEdgei].first();
150
151 tuple.procNo(Pstream::myProcNo());
152 tuple.faPatchi(patchId); // Tag as finiteArea patch
153 tuple.patchEdgei(patchEdgei);
154 tuple.meshFacei(meshFacei);
155 }
156
157 // Temporary local bookkeeping
158 {
159 auto& pairing = patchPairings[bndEdgei];
160
161 pairing.clear(); // invalidate
162 pairing.insert(patchId); // 'hash' into first location
163 pairing.patchEdgei_ = patchEdgei;
164 pairing.meshFacei_ = meshFacei;
165 }
166 }
167
168 if ((nBadEdges = returnReduce(badEdges.size(), sumOp<label>())) > 0)
169 {
170 edgeList dumpEdges(patch().edges(), badEdges.sortedToc());
171
172 vtk::lineWriter writer
173 (
174 patch().localPoints(),
175 dumpEdges,
176 fileName
177 (
178 mesh().time().globalPath()
179 / ("faMesh-construct.nonManifoldEdges")
180 )
181 );
182
183 writer.writeGeometry();
184
185 // CellData
186 writer.beginCellData();
187 writer.writeProcIDs();
188
190 << "Boundary edges not singly connected: "
191 << nBadEdges << '/' << nBoundaryEdges << nl;
192
194 (
196 patch(),
197 badEdges.sortedToc()
198 );
199
201 << "(debug) wrote " << writer.output().name() << nl;
202
204 }
205 badEdges.clear();
206
207 // Pass 2:
208 // Add in first connecting neighbour patch for the boundary edges.
209 // Need to examine all possibly connecting (non-processor) neighbours,
210 // but only need to check their boundary edges.
211
212 label nMissing = patchPairings.size();
213
214 for (label patchi = 0; patchi < nNonProcessor; ++patchi)
215 {
216 if (!nMissing) break; // Early exit
217
218 const polyPatch& pp = pbm[patchi];
219
220 // Check boundary edges
221 for
222 (
223 label patchEdgei = pp.nInternalEdges();
224 patchEdgei < pp.nEdges();
225 ++patchEdgei
226 )
227 {
228 const label bndEdgei =
229 edgeToBoundaryIndex.lookup(pp.meshEdge(patchEdgei), -1);
230
231 if (bndEdgei != -1)
232 {
233 // Has a matching owner boundary edge
234
235 auto& pairing = patchPairings[bndEdgei];
236
237 // Add neighbour (patchId, patchEdgei, meshFacei)
238 // 'hash' into the second location
239 // which does not insert the same value twice
240 if (pairing.insert(patchi))
241 {
242 // The attached patch face. Should only be one!
243 const labelList& edgeFaces = pp.edgeFaces()[patchEdgei];
244
245 if (edgeFaces.size() != 1)
246 {
247 pairing.erase(patchi);
248 badEdges.insert(badEdges.size());
249 continue;
250 }
251
252 const label patchFacei = edgeFaces[0];
253 const label meshFacei = patchFacei + pp.start();
254
255 // The neighbour information
256 pairing.patchEdgei_ = patchEdgei;
257 pairing.meshFacei_ = meshFacei;
258
259 --nMissing;
260 if (!nMissing) break; // Early exit
261 }
262 }
263 }
264 }
265
266 if ((nBadEdges = returnReduce(badEdges.size(), sumOp<label>())) > 0)
267 {
269 << "Had " << nBadEdges
270 << " boundary edges with missing or multiple edge connections"
271 << abort(FatalError);
272 }
273
274 // Combine local bookkeeping into final list
275 badEdges.clear();
276 for (label bndEdgei = 0; bndEdgei < nBoundaryEdges; ++bndEdgei)
277 {
278 const auto& pairing = patchPairings[bndEdgei];
279 const label nbrPatchi = pairing.second();
280 const label nbrPatchEdgei = pairing.patchEdgei_;
281 const label nbrMeshFacei = pairing.meshFacei_;
282
283 if (nbrMeshFacei >= 0)
284 {
285 // Add into primary bookkeeping
286 auto& tuple = bndEdgeConnections[bndEdgei].second();
287
288 tuple.procNo(Pstream::myProcNo());
289 tuple.patchi(nbrPatchi);
290 tuple.patchEdgei(nbrPatchEdgei);
291 tuple.meshFacei(nbrMeshFacei);
292 }
293 else if (!Pstream::parRun())
294 {
295 badEdges.insert(nInternalEdges + bndEdgei);
296 }
297 }
298 }
299
300
301 // ~~~~~~
302 // Serial - can return already
303 // ~~~~~~
304 if (!Pstream::parRun())
305 {
306 // Verbose report of missing edges - in serial
307
308 nBadEdges = badEdges.size();
309 if (nBadEdges)
310 {
311 edgeList dumpEdges(patch().edges(), badEdges.sortedToc());
312
313 vtk::lineWriter writer
314 (
315 patch().localPoints(),
316 dumpEdges,
317 fileName
318 (
319 mesh().time().globalPath()
320 / ("faMesh-construct.invalidEdges")
321 )
322 );
323
324 writer.writeGeometry();
325
326 // CellData
327 writer.beginCellData();
328 writer.writeProcIDs();
329
331 << "Boundary edges with missing/invalid neighbours: "
332 << nBadEdges << '/' << nBoundaryEdges << nl;
333
335 (
337 patch(),
338 badEdges.sortedToc()
339 );
340
342 << "(debug) wrote " << writer.output().name() << nl;
343
345 }
346
347 // Globally consistent ordering
348 patchTuple::sort(bndEdgeConnections);
349 return bndEdgeConnections;
350 }
351
352
353 // ~~~~~~~~
354 // Parallel
355 // ~~~~~~~~
356
358 << "Creating global coupling data" << endl;
359
360 const globalMeshData& globalData = mesh().globalData();
361 const indirectPrimitivePatch& cpp = globalData.coupledPatch();
362 const mapDistribute& map = globalData.globalEdgeSlavesMap();
363 const label nCoupledEdges = cpp.nEdges();
364
365 // Construct coupled edge usage with all data
366 List<bool> coupledEdgesUsed(map.constructSize(), false);
367
368 // Markup finiteArea boundary edges that are coupled across processors
369 for (label cppEdgei = 0; cppEdgei < nCoupledEdges; ++cppEdgei)
370 {
371 coupledEdgesUsed[cppEdgei] =
372 edgeToBoundaryIndex.found(cpp.meshEdge(cppEdgei));
373 }
374
376 << "Starting sync of boundary edge topology" << endl;
377
379 (
380 coupledEdgesUsed,
381 globalData.globalEdgeSlaves(),
382 globalData.globalEdgeTransformedSlaves(), // probably not used
383 map,
384 orEqOp<bool>()
385 );
386
387 if (debug)
388 {
389 label nAttached = 0;
390 for (label cppEdgei = 0; cppEdgei < nCoupledEdges; ++cppEdgei)
391 {
392 if (coupledEdgesUsed[cppEdgei])
393 {
394 ++nAttached;
395 }
396 }
397
399 << "Approx "
400 << returnReduce(nAttached, sumOp<label>())
401 << " connected boundary edges (total, some duplicates)" << endl;
402 }
403
404 // Combine information
405
406 // Normally 0 or 2 connections
407 List<DynamicList<patchTuple, 2>> gatheredConnections(map.constructSize());
408
409 // Map edges (mesh numbering) back to a coupled index in use
410 EdgeMap<label> edgeToCoupledIndex(2*nCoupledEdges);
411
412 // Pass 1
413 // Look for attached boundary edges
414 // - boundary edges from this side go into gathered connections
415 // - boundary edges connected from the other side are noted for later
416
417 for (label cppEdgei = 0; cppEdgei < nCoupledEdges; ++cppEdgei)
418 {
419 if (coupledEdgesUsed[cppEdgei])
420 {
421 const edge meshEdge(cpp.meshEdge(cppEdgei));
422
423 const label bndEdgei =
424 edgeToBoundaryIndex.lookup(meshEdge, -1);
425
426 if (bndEdgei != -1)
427 {
428 // A boundary finiteEdge edge (known from this side)
429
430 auto& gathered = gatheredConnections[cppEdgei];
431 gathered.setCapacity(2);
432 gathered.resize(1);
433 auto& tuple = gathered.last();
434
435 tuple = bndEdgeConnections[bndEdgei].first();
436 }
437 else
438 {
439 // Boundary edge connected from the other side
440 // - mark for it to be added later
441
442 edgeToCoupledIndex.insert(meshEdge, cppEdgei);
443 }
444 }
445 }
446
447 // Pass 2
448 // - add previously noted boundary edges (connected from other side)
449 // into gathered connections
450
451 badEdges.clear();
452 for (label patchi = 0; patchi < nNonProcessor; ++patchi)
453 {
454 if (edgeToCoupledIndex.empty()) break; // Early exit
455
456 const polyPatch& pp = pbm[patchi];
457
458 // Check boundary edges
459 for
460 (
461 label patchEdgei = pp.nInternalEdges();
462 patchEdgei < pp.nEdges();
463 ++patchEdgei
464 )
465 {
466 const edge meshEdge(pp.meshEdge(patchEdgei));
467
468 const label cppEdgei =
469 edgeToCoupledIndex.lookup(meshEdge, -1);
470
471 if (cppEdgei != -1)
472 {
473 // A known connection
474
475 // The attached patch face. Should only be one!
476 const labelList& edgeFaces = pp.edgeFaces()[patchEdgei];
477
478 if (edgeFaces.size() != 1)
479 {
480 badEdges.insert(cppEdgei);
481 continue;
482 }
483
484 const label patchFacei = edgeFaces[0];
485 const label meshFacei = patchFacei + pp.start();
486
487 auto& gathered = gatheredConnections[cppEdgei];
488 gathered.setCapacity(2);
489 gathered.resize(1);
490 auto& tuple = gathered.last();
491
492 tuple.procNo(Pstream::myProcNo());
493 tuple.patchi(patchi);
494 tuple.patchEdgei(patchEdgei);
495 tuple.meshFacei(meshFacei);
496
497 // Do not consider again
498 edgeToCoupledIndex.erase(meshEdge);
499
500 if (edgeToCoupledIndex.empty()) break; // Early exit
501 }
502 }
503 }
504
505 if ((nBadEdges = returnReduce(badEdges.size(), sumOp<label>())) > 0)
506 {
508 << "Had " << nBadEdges << " coupled boundary edges"
509 << " with missing or multiple edge connections"
510 << abort(FatalError);
511 }
512
514 << "Starting sync of boundary edge information" << endl;
515
517 (
518 gatheredConnections,
519 globalData.globalEdgeSlaves(),
520 globalData.globalEdgeTransformedSlaves(), // probably not used
521 map,
522 ListOps::appendEqOp<patchTuple>()
523 );
524
525
527 << "Collating sync information" << endl;
528
529 // Pick out gathered connections and add into primary bookkeeping
530 for (label cppEdgei = 0; cppEdgei < nCoupledEdges; ++cppEdgei)
531 {
532 const auto& gathered = gatheredConnections[cppEdgei];
533
534 const label bndEdgei =
535 edgeToBoundaryIndex.lookup(cpp.meshEdge(cppEdgei), -1);
536
537 if
538 (
539 // A boundary finiteEdge edge (known from this side)
540 bndEdgei != -1
541
542 // Gathered a one-to-one connection
543 && gathered.size() == 2
544 )
545 {
546 const auto& a = gathered[0];
547 const auto& b = gathered[1];
548
549 // Copy second side of connection
550 auto& connection = bndEdgeConnections[bndEdgei];
551
552 connection.second() = (connection.first() == b) ? a : b;
553 }
554 }
555
556 // Check missing/invalid
557 badEdges.clear();
558 for (label bndEdgei = 0; bndEdgei < nBoundaryEdges; ++bndEdgei)
559 {
560 const auto& connection = bndEdgeConnections[bndEdgei];
561
562 if (!connection.second().valid())
563 {
564 badEdges.insert(nInternalEdges + bndEdgei);
565 }
566 }
567
568
569 if (debug & 8)
570 {
571 // Boundary edges
572 {
573 vtk::lineWriter writer
574 (
575 patch().localPoints(),
576 patch().boundaryEdges(),
577 fileName
578 (
579 mesh().time().globalPath()
580 / ("faMesh-construct.boundaryEdges")
581 )
582 );
583
584 writer.writeGeometry();
585
586 // CellData
587 writer.beginCellData();
588 writer.writeProcIDs();
589
590 // For each boundary edge - the associate neighbour patch
591 labelList neighProc(nBoundaryEdges);
592 labelList neighPatch(nBoundaryEdges);
593 for (label bndEdgei = 0; bndEdgei < nBoundaryEdges; ++bndEdgei)
594 {
595 const auto& connection = bndEdgeConnections[bndEdgei];
596
597 neighProc[bndEdgei] = connection.second().procNo();
598 neighPatch[bndEdgei] = connection.second().patchi();
599 }
600
601 writer.write("neighProc", neighProc);
602 writer.write("neighPatch", neighPatch);
603 }
604
605 // finiteArea
606 {
608 (
609 patch(),
610 fileName
611 (
612 mesh().time().globalPath()
613 / ("faMesh-construct.faPatch")
614 )
615 );
616
617 writer.writeGeometry();
618
619 // CellData
620 writer.beginCellData();
621 writer.writeProcIDs();
622 }
623 }
624
625 // Verbose report of missing edges
626 if ((nBadEdges = returnReduce(badEdges.size(), sumOp<label>())) > 0)
627 {
628 edgeList dumpEdges(patch().edges(), badEdges.sortedToc());
629
630 vtk::lineWriter writer
631 (
632 patch().localPoints(),
633 dumpEdges,
634 fileName
635 (
636 mesh().time().globalPath()
637 / ("faMesh-construct.invalidEdges")
638 )
639 );
640
641 writer.writeGeometry();
642
643 // CellData
644 writer.beginCellData();
645 writer.writeProcIDs();
646
648 << "Boundary edges with missing/invalid neighbours: "
649 << nBadEdges << '/' << nBoundaryEdges << nl;
650
652 (
654 patch(),
655 badEdges.sortedToc()
656 );
657
659 << "(debug) wrote " << writer.output().name() << nl;
660
662 }
663
664
665 // Globally consistent ordering
666 patchTuple::sort(bndEdgeConnections);
667
669 << "Return sorted list of boundary connections" << endl;
670
671 return bndEdgeConnections;
672}
673
674
675void Foam::faMesh::setBoundaryConnections
676(
677 const List<Pair<patchTuple>>& bndEdgeConnections
678) const
679{
680 const label nInternalEdges = patch().nInternalEdges();
681 const label nBoundaryEdges = patch().nBoundaryEdges();
682
683 if (bndEdgeConnections.size() != nBoundaryEdges)
684 {
686 << "Sizing mismatch. Expected " << nBoundaryEdges
687 << " boundary edge connections, but had "
688 << bndEdgeConnections.size() << nl
689 << abort(FatalError);
690 }
691
692 bndConnectPtr_.reset
693 (
694 new List<labelPair>(nBoundaryEdges, labelPair(-1,-1))
695 );
696 auto& bndConnect = *bndConnectPtr_;
697
698 for (const auto& connection : bndEdgeConnections)
699 {
700 const auto& a = connection.first();
701 const auto& b = connection.second();
702
703 if (a.is_finiteArea() && a.is_localProc())
704 {
705 const label bndEdgei = (a.patchEdgei() - nInternalEdges);
706
707 bndConnect[bndEdgei].first() = b.procNo();
708 bndConnect[bndEdgei].second() = b.meshFacei();
709 }
710 else if (b.is_finiteArea() && b.is_localProc())
711 {
712 const label bndEdgei = (b.patchEdgei() - nInternalEdges);
713
714 bndConnect[bndEdgei].first() = a.procNo();
715 bndConnect[bndEdgei].second() = a.meshFacei();
716 }
717 else
718 {
720 << "Unexpected pairing input " << connection
721 << " ... programming error" << nl
722 << abort(FatalError);
723 }
724 }
725
726 label nInvalid = 0;
727 for (const auto& connection : bndConnect)
728 {
729 if (connection.first() < 0 || connection.second() < 0)
730 {
731 ++nInvalid;
732 }
733 }
734
735 if (Pstream::parRun())
736 {
737 reduce(nInvalid, sumOp<label>());
738 }
739
740 if (nInvalid)
741 {
743 << "Did not properly match " << nInvalid
744 << " boundary edges" << nl
745 << abort(FatalError);
746 }
747}
748
749
750void Foam::faMesh::calcBoundaryConnections() const
751{
752 setBoundaryConnections(this->getBoundaryEdgeConnections());
753}
754
755
756// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
757
759{
760 const auto& connections = this->boundaryConnections();
761
762 labelHashSet procsUsed(2*Pstream::nProcs());
763
764 for (const labelPair& tuple : connections)
765 {
766 procsUsed.insert(tuple.first());
767 }
768
769 procsUsed.erase(-1); // placeholder value
770 procsUsed.erase(Pstream::myProcNo());
771
772 return procsUsed.sortedToc();
773}
774
775
777{
778 const auto& connections = this->boundaryConnections();
779
780 Map<label> procCount(2*Pstream::nProcs());
781
782 for (const labelPair& tuple : connections)
783 {
784 ++procCount(tuple.first());
785 }
786 procCount.erase(-1); // placeholder value
787 procCount.erase(Pstream::myProcNo());
788
789 // Flatten as list
790 List<labelPair> output(procCount.size());
791 label count = 0;
792 for (const label proci : procCount.sortedToc())
793 {
794 output[count].first() = proci;
795 output[count].second() = procCount[proci]; // size
796 ++count;
797 }
798
799 return output;
800}
801
802
804{
805 if (!haloMapPtr_)
806 {
807 haloMapPtr_.reset(new faMeshBoundaryHalo(*this));
808 }
809
810 return *haloMapPtr_;
811}
812
813
814void Foam::faMesh::calcHaloFaceGeometry() const
815{
816 if (haloFaceCentresPtr_ || haloFaceNormalsPtr_)
817 {
819 << "Halo centres/normals already calculated"
820 << exit(FatalError);
821 }
822
824 << "Calculating halo face centres/normals" << endl;
825
826 const faceList& faces = mesh().faces();
827 const pointField& points = mesh().points();
828
829 const faMeshBoundaryHalo& halo = boundaryHaloMap();
830
831 const labelList& inputFaceIds = halo.inputMeshFaces();
832
833 haloFaceCentresPtr_.reset(new pointField);
834 haloFaceNormalsPtr_.reset(new vectorField);
835
836 auto& centres = *haloFaceCentresPtr_;
837 auto& normals = *haloFaceNormalsPtr_;
838
839 centres.resize(inputFaceIds.size());
840 normals.resize(inputFaceIds.size());
841
842 // My values
843 forAll(inputFaceIds, i)
844 {
845 const face& f = faces[inputFaceIds[i]];
846
847 centres[i] = f.centre(points);
848 normals[i] = f.unitNormal(points);
849 }
850
851 // Swap information and resize
852 halo.distributeSparse(centres);
853 halo.distributeSparse(normals);
854}
855
856
858{
859 if (!haloFaceCentresPtr_ || !haloFaceNormalsPtr_)
860 {
861 calcHaloFaceGeometry();
862 }
863
864 return *haloFaceCentresPtr_;
865}
866
867
869{
870 if (!haloFaceCentresPtr_ || !haloFaceNormalsPtr_)
871 {
872 calcHaloFaceGeometry();
873 }
874
875 return *haloFaceNormalsPtr_;
876}
877
878
880Foam::faMesh::haloFaceCentres(const label patchi) const
881{
882 if (patchi < 0 || patchi >= boundary().size())
883 {
885 << "Patch " << patchi << " is out-of-range 0.."
886 << (boundary().size()-1) << nl
887 << exit(FatalError);
888 }
889
891 (
893 (
894 boundarySubset
895 (
896 this->haloFaceCentres(),
897 boundary()[patchi].edgeLabels()
898 )
899 )
900 );
901}
902
903
905Foam::faMesh::haloFaceNormals(const label patchi) const
906{
907 if (patchi < 0 || patchi >= boundary().size())
908 {
910 << "Patch " << patchi << " is out-of-range 0.."
911 << (boundary().size()-1) << nl
912 << exit(FatalError);
913 }
914
916 (
918 (
919 boundarySubset
920 (
921 this->haloFaceNormals(),
922 boundary()[patchi].edgeLabels()
923 )
924 )
925 );
926}
927
928
929// ************************************************************************* //
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:440
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
label nBoundaryEdges() const
Number of boundary edges == (nEdges() - nInternalEdges())
label nInternalEdges() const
Number of internal edges.
edge meshEdge(const label edgei) const
From patch edge to global edge using meshPoints.
const labelListList & edgeFaces() const
Return edge-face addressing.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
static int debug
Debug switch.
Definition: data.H:77
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
void clear()
'Clears' edge by setting both ends to invalid point labels.
Definition: edgeI.H:231
void sort()
Sort element lists numerically.
Definition: ensightCells.C:149
Class for obtaining halo face data for the boundary edges. The ordering follows that natural edge ord...
void reset(const faMesh &mesh)
Redefine map and connectivity for a mesh.
const vectorField & haloFaceNormals() const
Face normals of boundary halo neighbours.
const Time & time() const
Return reference to time.
Definition: faMesh.C:673
const pointField & haloFaceCentres() const
Face centres of boundary halo neighbours.
label nBoundaryEdges() const noexcept
Number of boundary edges (== nEdges - nInternalEdges)
Definition: faMeshI.H:74
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMeshI.H:31
const edgeList & edges() const noexcept
Return local edges with reordered boundary.
Definition: faMeshI.H:92
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMeshI.H:128
label nInternalEdges() const noexcept
Number of internal faces.
Definition: faMeshI.H:68
const faGlobalMeshData & globalData() const
Return parallel info.
Definition: faMesh.C:893
labelList boundaryProcs() const
List< labelPair > boundaryProcSizes() const
const faMeshBoundaryHalo & boundaryHaloMap() const
Mapping/swapping for boundary halo neighbours.
static void syncData(List< Type > &elems, const labelListList &slaves, const labelListList &transformedSlaves, const mapDistribute &slavesMap, const globalIndexAndTransform &, const CombineOp &cop, const TransformOp &top)
Helper: synchronise data with transforms.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
label nInternalFaces() const noexcept
Number of internal faces.
int myProcNo() const noexcept
Return processor number.
A class for managing temporary objects.
Definition: tmp.H:65
volScalarField & p
patchWriters clear()
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const label nNonProcessor
const pointField & points
label patchId(-1)
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
const std::string patch
OpenFOAM patch number as a std::string.
GenericPatchWriter< uindirectPrimitivePatch > uindirectPatchWriter
Write uindirectPrimitivePatch faces/points (optionally with fields) as a vtp file or a legacy vtk fil...
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:57
List< label > labelList
A List of labels.
Definition: List.H:66
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
static void printPatchEdges(Ostream &os, const PatchType &p, const labelList &edgeIds, label maxOutput=10)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Field< vector > vectorField
Specialisation of Field<T> for vector.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
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 f(nPoints)
volScalarField & b
Definition: createFields.H:27
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333