faBoundaryMesh.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) 2018-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 "faBoundaryMesh.H"
30#include "faMesh.H"
31#include "globalIndex.H"
32#include "primitiveMesh.H"
33#include "processorFaPatch.H"
34#include "PtrListOps.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
46bool Foam::faBoundaryMesh::hasGroupIDs() const
47{
48 if (groupIDsPtr_)
49 {
50 // Use existing cache
51 return !groupIDsPtr_->empty();
52 }
53
54 const faPatchList& patches = *this;
55
56 for (const faPatch& p : patches)
57 {
58 if (!p.inGroups().empty())
59 {
60 return true;
61 }
62 }
63
64 return false;
65}
66
67
68void Foam::faBoundaryMesh::calcGroupIDs() const
69{
70 if (groupIDsPtr_)
71 {
72 return; // Or FatalError
73 }
74
75 groupIDsPtr_.reset(new HashTable<labelList>(16));
76 auto& groupLookup = *groupIDsPtr_;
77
78 const faPatchList& patches = *this;
79
80 forAll(patches, patchi)
81 {
82 const wordList& groups = patches[patchi].inGroups();
83
84 for (const word& groupName : groups)
85 {
86 groupLookup(groupName).append(patchi);
87 }
88 }
89
90 // Remove groups that clash with patch names
91 forAll(patches, patchi)
92 {
93 if (groupLookup.erase(patches[patchi].name()))
94 {
96 << "Removed group '" << patches[patchi].name()
97 << "' which clashes with patch " << patchi
98 << " of the same name."
99 << endl;
100 }
101 }
102}
103
104
105bool Foam::faBoundaryMesh::readContents(const bool allowReadIfPresent)
106{
107 if
108 (
109 readOpt() == IOobject::MUST_READ
110 || readOpt() == IOobject::MUST_READ_IF_MODIFIED
111 ||
112 (
113 allowReadIfPresent
114 && (readOpt() == IOobject::READ_IF_PRESENT && headerOk())
115 )
116 )
117 {
118 // Warn for MUST_READ_IF_MODIFIED
119 warnNoRereading<faBoundaryMesh>();
120
121 faPatchList& patches = *this;
122
123 // Read faPatch list
124 Istream& is = readStream(typeName);
125
126 // Read patches as entries
127 PtrList<entry> patchEntries(is);
128 patches.resize(patchEntries.size());
129
130 // Transcribe
131 forAll(patches, patchi)
132 {
134 (
135 patchi,
137 (
138 patchEntries[patchi].keyword(),
139 patchEntries[patchi].dict(),
140 patchi,
141 *this
142 )
143 );
144 }
145
146 is.check(FUNCTION_NAME);
147 close();
148 return true;
149 }
150
151 return false;
152}
153
154
155// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
156
158(
159 const IOobject& io,
160 const faMesh& mesh
161)
162:
163 faPatchList(),
165 mesh_(mesh)
166{
167 readContents(false); // READ_IF_PRESENT allowed: False
168}
169
170
172(
173 const IOobject& io,
174 const faMesh& pm,
175 const label size
176)
177:
178 faPatchList(size),
180 mesh_(pm)
181{}
182
183
184// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185
187{
188 // processorFaPatch geometry triggers calculation of pointNormals.
189 // This uses parallel comms and hence will not be trigggered
190 // on processors that do not have a processorFaPatch so instead
191 // force construction.
192 (void)mesh_.pointAreaNormals();
193
195
196 if
197 (
200 )
201 {
202 forAll(*this, patchi)
203 {
204 operator[](patchi).initGeometry(pBufs);
205 }
206
207 pBufs.finishedSends();
208
209 forAll(*this, patchi)
210 {
211 operator[](patchi).calcGeometry(pBufs);
212 }
213 }
214 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
215 {
216 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
217
218 // Dummy.
219 pBufs.finishedSends();
220
221 for (const auto& patchEval : patchSchedule)
222 {
223 const label patchi = patchEval.patch;
224
225 if (patchEval.init)
226 {
227 operator[](patchi).initGeometry(pBufs);
228 }
229 else
230 {
231 operator[](patchi).calcGeometry(pBufs);
232 }
233 }
234 }
235}
236
237
240{
241 const faPatchList& patches = *this;
242
244
245 forAll(list, patchi)
246 {
247 list.set(patchi, &patches[patchi].edgeLabels());
248 }
249
250 return list;
251}
252
253
256{
257 const faPatchList& patches = *this;
258
260
261 forAll(list, patchi)
262 {
263 list.set(patchi, &patches[patchi].edgeFaces());
264 }
265
266 return list;
267}
268
269
271{
272 const faPatchList& patches = *this;
273
275
276 forAll(list, patchi)
277 {
278 const lduInterface* lduPtr = isA<lduInterface>(patches[patchi]);
279
280 if (lduPtr)
281 {
282 list.set(patchi, lduPtr);
283 }
284 }
285
286 return list;
287}
288
289
291{
292 const faPatchList& patches = *this;
293
294 label nonProc = 0;
295
296 for (const faPatch& p : patches)
297 {
298 if (isA<processorFaPatch>(p))
299 {
300 break;
301 }
302
303 ++nonProc;
304 }
305
306 return nonProc;
307}
308
309
312{
313 if (!groupIDsPtr_)
314 {
315 calcGroupIDs();
316 }
317
318 return *groupIDsPtr_;
319}
320
321
323(
324 const word& groupName,
325 const labelUList& patchIDs
326)
327{
328 groupIDsPtr_.clear();
329
330 faPatchList& patches = *this;
331
332 boolList donePatch(patches.size(), false);
333
334 // Add to specified patches
335 for (const label patchi : patchIDs)
336 {
337 patches[patchi].inGroups().appendUniq(groupName);
338 donePatch[patchi] = true;
339 }
340
341 // Remove from other patches
342 forAll(patches, patchi)
343 {
344 if (!donePatch[patchi])
345 {
346 wordList& groups = patches[patchi].inGroups();
347
348 if (groups.found(groupName))
349 {
350 label newi = 0;
351 forAll(groups, i)
352 {
353 if (groups[i] != groupName)
354 {
355 groups[newi++] = groups[i];
356 }
357 }
358 groups.resize(newi);
359 }
360 }
361 }
362}
363
364
366{
367 return PtrListOps::get<word>(*this, nameOp<faPatch>());
368}
369
370
372{
373 return PtrListOps::get<word>(*this, typeOp<faPatch>());
374}
375
376
378{
379 // Manually: faPatch does not have independent start() information
380
381 const faPatchList& patches = *this;
382
383 labelList list(patches.size());
384
385 label beg = mesh_.nInternalEdges();
386 forAll(patches, patchi)
387 {
388 const label len = patches[patchi].nEdges();
389 list[patchi] = beg;
390 beg += len;
391 }
392 return list;
393}
394
395
397{
398 return
399 PtrListOps::get<label>
400 (
401 *this,
402 [](const faPatch& p) { return p.nEdges(); } // avoid virtual
403 );
404}
405
406
408{
409 const faPatchList& patches = *this;
410
412
413 label beg = mesh_.nInternalEdges();
414 forAll(patches, patchi)
415 {
416 const label len = patches[patchi].nEdges();
417 list[patchi].reset(beg, len);
418 beg += len;
419 }
420 return list;
421}
422
423
425{
426 return mesh_.nInternalEdges();
427}
428
429
431{
432 return mesh_.nBoundaryEdges();
433}
434
435
437{
438 return labelRange(mesh_.nInternalEdges(), mesh_.nBoundaryEdges());
439}
440
441
443(
444 const wordRe& matcher,
445 const bool useGroups
446) const
447{
448 if (matcher.empty())
449 {
450 return labelList();
451 }
452
453 // Only check groups if requested and they exist
454 const bool checkGroups = (useGroups && this->hasGroupIDs());
455
456 labelHashSet ids;
457
458 if (matcher.isPattern())
459 {
460 if (checkGroups)
461 {
462 const auto& groupLookup = groupPatchIDs();
463 forAllConstIters(groupLookup, iter)
464 {
465 if (matcher.match(iter.key()))
466 {
467 // Hash ids associated with the group
468 ids.insert(iter.val());
469 }
470 }
471 }
472
473 if (ids.empty())
474 {
475 return PtrListOps::findMatching(*this, matcher);
476 }
477 else
478 {
479 ids.insert(PtrListOps::findMatching(*this, matcher));
480 }
481 }
482 else
483 {
484 // Literal string.
485 // Special version of above for reduced memory footprint
486
487 const label patchId = PtrListOps::firstMatching(*this, matcher);
488
489 if (patchId >= 0)
490 {
491 return labelList(one{}, patchId);
492 }
493 else if (checkGroups)
494 {
495 const auto iter = groupPatchIDs().cfind(matcher);
496
497 if (iter.found())
498 {
499 // Hash ids associated with the group
500 ids.insert(iter.val());
501 }
502 }
503 }
504
505 return ids.sortedToc();
506}
507
508
510(
511 const wordRes& matcher,
512 const bool useGroups
513) const
514{
515 if (matcher.empty())
516 {
517 return labelList();
518 }
519 else if (matcher.size() == 1)
520 {
521 return this->indices(matcher.first(), useGroups);
522 }
523
524 labelHashSet ids;
525
526 // Only check groups if requested and they exist
527 if (useGroups && this->hasGroupIDs())
528 {
529 ids.resize(2*this->size());
530
531 const auto& groupLookup = groupPatchIDs();
532 forAllConstIters(groupLookup, iter)
533 {
534 if (matcher.match(iter.key()))
535 {
536 // Hash ids associated with the group
537 ids.insert(iter.val());
538 }
539 }
540 }
541
542 if (ids.empty())
543 {
544 return PtrListOps::findMatching(*this, matcher);
545 }
546 else
547 {
548 ids.insert(PtrListOps::findMatching(*this, matcher));
549 }
550
551 return ids.sortedToc();
552}
553
554
555Foam::label Foam::faBoundaryMesh::findIndex(const wordRe& key) const
556{
557 if (key.empty())
558 {
559 return -1;
560 }
561 return PtrListOps::firstMatching(*this, key);
562}
563
564
566(
567 const word& patchName,
568 bool allowNotFound
569) const
570{
571 if (patchName.empty())
572 {
573 return -1;
574 }
575
576 const label patchId = PtrListOps::firstMatching(*this, patchName);
577
578 if (patchId >= 0)
579 {
580 return patchId;
581 }
582
583 if (!allowNotFound)
584 {
586 << "Patch '" << patchName << "' not found. "
587 << "Available patch names: " << names() << endl
588 << exit(FatalError);
589 }
590
591 // Patch not found
592 if (debug)
593 {
594 Pout<< "label faBoundaryMesh::findPatchID(const word&) const"
595 << "Patch named " << patchName << " not found. "
596 << "Available patch names: " << names() << endl;
597 }
598
599 // Not found, return -1
600 return -1;
601}
602
603
604Foam::label Foam::faBoundaryMesh::whichPatch(const label edgeIndex) const
605{
606 // Find out which patch the current face belongs to by comparing label
607 // with patch start labels.
608 // If the face is internal, return -1;
609 // if it is off the end of the list, abort
610 if (edgeIndex < mesh().nInternalEdges())
611 {
612 return -1;
613 }
614 else if (edgeIndex >= mesh().nEdges())
615 {
617 << "Edge " << edgeIndex
618 << " out of bounds. Number of geometric edges " << mesh().nEdges()
619 << abort(FatalError);
620 }
621
622 forAll(*this, patchi)
623 {
624 label start = mesh_.patchStarts()[patchi];
625 label size = operator[](patchi).faPatch::size();
626
627 if
628 (
629 edgeIndex >= start
630 && edgeIndex < start + size
631 )
632 {
633 return patchi;
634 }
635 }
636
637 // If not in any of above, it's trouble!
639 << "error in patch search algorithm"
640 << abort(FatalError);
641
642 return -1;
643}
644
645
646bool Foam::faBoundaryMesh::checkParallelSync(const bool report) const
647{
648 if (!Pstream::parRun())
649 {
650 return false;
651 }
652
653 const faBoundaryMesh& bm = *this;
654
655 bool hasError = false;
656
657 // Collect non-proc patches and check proc patches are last.
658 wordList localNames(bm.size());
659 wordList localTypes(bm.size());
660
661 label nonProci = 0;
662
663 forAll(bm, patchi)
664 {
665 if (!isA<processorFaPatch>(bm[patchi]))
666 {
667 if (nonProci != patchi)
668 {
669 // A processor patch in between normal patches!
670 hasError = true;
671
672 if (debug || report)
673 {
674 Pout<< " ***Problem with boundary patch " << patchi
675 << " name:" << bm[patchi].name()
676 << " type:" << bm[patchi].type()
677 << " - seems to be preceeded by processor patches."
678 << " This is usually a problem." << endl;
679 }
680 }
681 else
682 {
683 localNames[nonProci] = bm[patchi].name();
684 localTypes[nonProci] = bm[patchi].type();
685 ++nonProci;
686 }
687 }
688 }
689 localNames.resize(nonProci);
690 localTypes.resize(nonProci);
691
692 // Check and report error(s) on master
693
695 (
696 // Don't need to collect master itself
697 (Pstream::master() ? 0 : nonProci),
699 );
700
701 const wordList allNames(procAddr.gather(localNames));
702 const wordList allTypes(procAddr.gather(localTypes));
703
704 // Automatically restricted to master
705 for (const int proci : procAddr.subProcs())
706 {
707 const auto procNames(allNames.slice(procAddr.range(proci)));
708 const auto procTypes(allTypes.slice(procAddr.range(proci)));
709
710 if (procNames != localNames || procTypes != localTypes)
711 {
712 hasError = true;
713
714 if (debug || report)
715 {
716 Info<< " ***Inconsistent patches across processors, "
717 "processor0 has patch names:" << localNames
718 << " patch types:" << localTypes
719 << " processor" << proci
720 << " has patch names:" << procNames
721 << " patch types:" << procTypes
722 << endl;
723 }
724 }
725 }
726
727 // Reduce (not broadcast) to respect local out-of-order errors (first loop)
728 reduce(hasError, orOp<bool>());
729
730 return hasError;
731}
732
733
734bool Foam::faBoundaryMesh::checkDefinition(const bool report) const
735{
736 label nextPatchStart = mesh().nInternalEdges();
737 const faBoundaryMesh& bm = *this;
738
739 bool hasError = false;
740
741 forAll(bm, patchi)
742 {
743 if (bm[patchi].start() != nextPatchStart && !hasError)
744 {
745 hasError = true;
746
748 << " ****Problem with boundary patch " << patchi
749 << " named " << bm[patchi].name()
750 << " of type " << bm[patchi].type()
751 << ". The patch should start on face no " << nextPatchStart
752 << " and the patch specifies " << bm[patchi].start()
753 << "." << endl
754 << "Possibly consecutive patches have this same problem."
755 << " Suppressing future warnings." << endl;
756 }
757
758 // Warn about duplicate boundary patches?
759
760 nextPatchStart += bm[patchi].faPatch::size();
761 }
762
763 if (hasError)
764 {
766 << "This mesh is not valid: boundary definition is in error."
767 << endl;
768 }
769 else
770 {
771 if (debug || report)
772 {
773 Info << "Boundary definition OK." << endl;
774 }
775 }
776
777 return hasError;
778}
779
780
782{
783 // processorFaPatch geometry triggers calculation of pointNormals.
784 // This uses parallel comms and hence will not be trigggered
785 // on processors that do not have a processorFaPatch so instead
786 // force construction.
787 (void)mesh_.pointAreaNormals();
788
790
791 if
792 (
795 )
796 {
797 forAll(*this, patchi)
798 {
799 operator[](patchi).initMovePoints(pBufs, p);
800 }
801
802 pBufs.finishedSends();
803
804 forAll(*this, patchi)
805 {
806 operator[](patchi).movePoints(pBufs, p);
807 }
808 }
809 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
810 {
811 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
812
813 // Dummy.
814 pBufs.finishedSends();
815
816 for (const auto& schedEval : patchSchedule)
817 {
818 const label patchi = schedEval.patch;
819
820 if (schedEval.init)
821 {
822 operator[](patchi).initMovePoints(pBufs, p);
823 }
824 else
825 {
826 operator[](patchi).movePoints(pBufs, p);
827 }
828 }
829 }
830}
831
832
834{
836
837 if
838 (
841 )
842 {
843 forAll(*this, patchi)
844 {
845 operator[](patchi).initUpdateMesh(pBufs);
846 }
847
848 pBufs.finishedSends();
849
850 forAll(*this, patchi)
851 {
852 operator[](patchi).updateMesh(pBufs);
853 }
854 }
855 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
856 {
857 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
858
859 // Dummy.
860 pBufs.finishedSends();
861
862 for (const auto& schedEval : patchSchedule)
863 {
864 const label patchi = schedEval.patch;
865
866 if (schedEval.init)
867 {
868 operator[](patchi).initUpdateMesh(pBufs);
869 }
870 else
871 {
872 operator[](patchi).updateMesh(pBufs);
873 }
874 }
875 }
876}
877
878
880{
881 const faPatchList& patches = *this;
882
884
885 for (const faPatch& p : patches)
886 {
887 os.beginBlock(p.name());
888 os << p;
889 os.endBlock();
890 }
891
893
895 return os.good();
896}
897
898
900(
901 IOstreamOption streamOpt,
902 const bool valid
903) const
904{
905 // Allow/disallow compression?
906 // 1. keep readable
907 // 2. save some space
908 // ??? streamOpt.compression(IOstreamOption::UNCOMPRESSED);
909 return regIOobject::writeObject(streamOpt, valid);
910}
911
912
913// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
914
916{
917 bm.writeData(os);
918 return os;
919}
920
921
922// ************************************************************************* //
Functions to operate on Pointer Lists.
globalIndex procAddr(aMesh.nFaces())
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTableI.H:59
void resize(const label sz)
Resize the hash table for efficiency.
Definition: HashTable.C:601
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:180
The IOstreamOption is a simple container for options an IOstream can normally have.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
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
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:105
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:87
Buffers for inter-processor communications streams (UOPstream, UIPstream).
UPstream::commsTypes commsType() const noexcept
The communications type of the stream.
void finishedSends(const bool wait=true)
Mark sends as done.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
void append(T *ptr)
Append an element to the end of the list.
Definition: PtrListI.H:113
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:103
const labelUList & indices() const noexcept
Return the list of sorted indices (updated every sort).
Definition: SortListI.H:54
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
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:265
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
Definition: SubList.H:165
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
const T * set(const label i) const
Definition: UPtrList.H:248
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Finite area boundary mesh.
UPtrList< const labelUList > edgeFaces() const
Return a list of edgeFaces for each patch.
bool checkDefinition(const bool report=false) const
Check boundary definition.
label nNonProcessor() const
The number of patches before the first processor patch.
void calcGeometry()
Calculate the geometry for the patches.
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
labelList patchStarts() const
Return a list of patch start indices.
lduInterfacePtrsList interfaces() const
labelRange range() const
The edge range for all boundary edges.
wordList types() const
Return a list of patch types.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
bool writeData(Ostream &os) const
The writeData member function required by regIOobject.
void setGroup(const word &groupName, const labelUList &patchIDs)
Set/add group with patches.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
List< labelRange > patchRanges() const
Return a list of patch ranges.
label whichPatch(const label edgeIndex) const
Return patch index for a given edge label.
labelList patchSizes() const
Return a list of patch sizes (number of edges in each patch)
wordList names() const
Return a list of patch names.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
void updateMesh()
Correct faBoundaryMesh after topology update.
UPtrList< const labelUList > edgeLabels() const
Return a list of edgeLabels for each patch.
label nEdges() const
The number of boundary edges for the underlying mesh.
label findIndex(const wordRe &key) const
Return patch index for the first match, return -1 if not found.
label start() const
The start label of the edges in the faMesh edges list.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:78
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
labelRange range() const
Return start/size range of local processor data.
Definition: globalIndexI.H:232
static void gather(const labelUList &offsets, const label comm, const ProcIDsContainer &procIDs, const UList< Type > &fld, List< Type > &allFld, const int tag=UPstream::msgType(), const UPstream::commsTypes=UPstream::commsTypes::nonBlocking)
Collect data in processor order on master (== procIDs[0]).
labelRange subProcs() const noexcept
Range of process indices for addressed sub-offsets (processes)
Definition: globalIndexI.H:159
const lduSchedule & patchSchedule() const noexcept
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:58
An abstract base class for implicitly-coupled interfaces e.g. processor and cyclic patches.
Definition: lduInterface.H:58
void movePoints()
Update for new mesh geometry.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
label nInternalEdges() const
Internal edges using 0,1 or 2 boundary points.
label nEdges() const
Number of mesh edges.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
splitCell * master() const
Definition: splitCell.H:113
@ BEGIN_LIST
Begin list [isseparator].
Definition: token.H:155
@ END_LIST
End list [isseparator].
Definition: token.H:156
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:83
bool isPattern() const noexcept
The wordRe is a pattern, not a literal string.
Definition: wordReI.H:111
bool match(const std::string &text, bool literal=false) const
Smart match as regular expression or as a string.
Definition: wordReI.H:200
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
bool match(const std::string &text, bool literal=false) const
Smart match as literal or regex, stopping on the first match.
Definition: wordResI.H:91
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
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
label patchId(-1)
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
#define InfoInFunction
Report an information message using Foam::Info.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
labelList findMatching(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Extract list indices for all items with 'name()' that matches.
label firstMatching(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Find first list item with 'name()' that matches, -1 on failure.
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:349
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)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition: faPatch.H:66
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:356
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:238
Extract type (as a word) from an object, typically using its type() method.
Definition: word.H:255