Loading...
Searching...
No Matches
pointBoundaryMesh.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2025 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 "pointBoundaryMesh.H"
30#include "polyBoundaryMesh.H"
31#include "facePointPatch.H"
32#include "pointMesh.H"
33#include "PstreamBuffers.H"
34#include "lduSchedule.H"
36#include "processorPointPatch.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
43}
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
47bool Foam::pointBoundaryMesh::hasGroupIDs() const
48{
49 if (groupIDsPtr_)
50 {
51 // Use existing cache
52 return !groupIDsPtr_->empty();
53 }
54
55 const auto& patches = *this;
56
57 for (const auto& p : patches)
58 {
59 if (!p.inGroups().empty())
60 {
61 return true;
62 }
63 }
64
65 return false;
66}
67
68
69void Foam::pointBoundaryMesh::calcGroupIDs() const
70{
71 if (groupIDsPtr_)
72 {
73 return; // Or FatalError
74 }
75
76 groupIDsPtr_.emplace(16);
77 auto& groupLookup = *groupIDsPtr_;
78
79 const auto& patches = *this;
80
81 forAll(patches, patchi)
82 {
83 for (const word& groupName : patches[patchi].inGroups())
84 {
85 groupLookup(groupName).push_back(patchi);
86 }
87 }
88
89 // Remove groups that clash with patch names
90 forAll(patches, patchi)
91 {
92 if (groupLookup.empty())
93 {
94 break; // Early termination
95 }
96 else if (groupLookup.erase(patches[patchi].name()))
97 {
99 << "Removed group '" << patches[patchi].name()
100 << "' which clashes with patch " << patchi
101 << " of the same name."
102 << endl;
103 }
104 }
105}
106
107
108void Foam::pointBoundaryMesh::addPatches(const polyBoundaryMesh& pbm)
109{
110 // Set boundary patches
111 pointPatchList& patches = *this;
112
114
115 forAll(pbm, patchi)
116 {
117 // NB: needs ptr() to get *pointPatch instead of *facePointPatch
118 patches.set(patchi, facePointPatch::New(pbm[patchi], *this).ptr());
119 }
120}
121
122
123// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124
125Foam::pointBoundaryMesh::pointBoundaryMesh
126(
127 const pointMesh& m,
128 const polyBoundaryMesh& pbm
129)
130:
132 regIOobject
133 (
134 IOobject
135 (
136 "boundary",
137 //m.thisDb().time().findInstance(m.meshDir(), "boundary"),
138 pbm.mesh().facesInstance(),
139 polyMesh::meshSubDir/pointMesh::meshSubDir,
140 m.thisDb(),
141 IOobject::NO_READ,
142 IOobject::NO_WRITE,
143 false // Avoid conflict with polyMesh/boundary
144 )
145 ),
146 mesh_(m)
147{
148 addPatches(pbm);
149
150 if (debug)
151 {
152 pointPatchList& Patches = *this;
153
154 Pout<< "pointBoundaryMesh::pointBoundaryMesh"
155 << "(const pointMesh&, const polyBoundaryMesh&): "
156 << "constructed pointBoundaryMesh:" << endl;
158 for (const auto& pp : Patches)
159 {
160 Pout<< indent
161 << "index:" << pp.index() << " patch:" << pp.name()
162 << " type:" << pp.type() << endl;
163 }
165 }
166}
167
168
169Foam::pointBoundaryMesh::pointBoundaryMesh
170(
171 const IOobject& io,
172 const pointMesh& m,
173 const polyBoundaryMesh& pbm
174)
175:
178 (
180 (
181 "boundary",
182 io.instance(),
183 polyMesh::meshSubDir/pointMesh::meshSubDir,
184 m.thisDb(),
185 io.readOpt(),
186 io.writeOpt(),
187 false //io.registerObject() // or always set to false?
188 )
189 ),
190 mesh_(m)
191{
192 pointPatchList& Patches = *this;
193
194 if (isReadRequired() || (isReadOptional() && headerOk()))
195 {
197 {
199 << "Specified IOobject::MUST_READ_IF_MODIFIED but class"
200 << " does not support automatic rereading."
201 << endl;
202 }
203
204 if (debug)
205 {
206 Pout<< "pointBoundaryMesh::pointBoundaryMesh"
207 << "(const IOobject&, const pointMesh&,"
208 << " const polyBoundaryMesh&): "
209 << "Constructing from boundary file " << objectRelPath()
210 << endl;
211 }
212
213 // Read pointPatchList
214 Istream& is = readStream(typeName);
215
216 PtrList<entry> patchEntries(is);
217 Patches.setSize(patchEntries.size());
218
219 forAll(Patches, patchi)
220 {
221 // Try construct-from-dictionary first
222 const word& name = patchEntries[patchi].keyword();
223
224 autoPtr<pointPatch> pPtr
225 (
227 (
228 name,
229 patchEntries[patchi].dict(),
230 patchi,
231 *this
232 )
233 );
234
235 if (!pPtr)
236 {
237 const label polyPatchi = pbm.findPatchID(name, false);
238 // Try as facePointPatch from polyPatch
239 pPtr = facePointPatch::New(pbm[polyPatchi], *this);
240 pPtr->index() = patchi;
241 }
242
243 Patches.set(patchi, pPtr);
244 }
245
246 // Check state of IOstream
247 is.check
248 (
249 "pointBoundaryMesh::pointBoundaryMesh"
250 "(const IOobject&, const pointMesh&,"
251 " const polyBoundaryMesh&)"
252 );
253
254 close();
255 }
256 else
257 {
258 if (debug)
259 {
260 Pout<< "pointBoundaryMesh::pointBoundaryMesh"
261 << "(const IOobject&, const pointMesh&,"
262 << " const polyBoundaryMesh&): "
263 << "Constructing from polyBoundaryMesh only"
264 << endl;
265 }
266
267 addPatches(pbm);
268 }
269
270
271 if (debug)
272 {
273 Pout<< "pointBoundaryMesh::pointBoundaryMesh"
274 << "(const IOobject&, const pointMesh&, const polyBoundaryMesh&): "
275 << "constructed pointBoundaryMesh:" << endl;
277 for (const auto& pp : Patches)
278 {
279 Pout<< indent
280 << "index:" << pp.index() << " patch:" << pp.name()
281 << " type:" << pp.type() << endl;
282 }
284 }
285}
286
287
288// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
289
291{
292 const pointPatchList& patches = *this;
293
294 label count = 0;
295
296 for (const auto& p : patches)
297 {
299 {
300 break;
301 }
302
304 }
305
306 return count;
307}
308
309
311{
312 const pointPatchList& patches = *this;
313
314 label count = 0;
315
316 for (const auto& p : patches)
317 {
319 {
320 ++count;
322 }
323
324 return count;
325}
326
331}
332
342 return
344 (
345 *this,
346 [](const pointPatch& p) { return p.physicalType(); }
347 );
348}
349
350
353{
354 if (!groupIDsPtr_)
355 {
356 calcGroupIDs();
357 }
358
359 return *groupIDsPtr_;
360}
361
362
364(
365 const wordRe& matcher,
366 const bool useGroups
367) const
368{
369 if (matcher.empty())
370 {
371 return labelList();
372 }
373
374 // Only check groups if requested and they exist
375 const bool checkGroups = (useGroups && this->hasGroupIDs());
376
377 labelHashSet ids;
378
379 if (matcher.isPattern())
380 {
381 if (checkGroups)
382 {
383 const auto& groupLookup = groupPatchIDs();
384 forAllConstIters(groupLookup, iter)
385 {
386 if (matcher(iter.key()))
387 {
388 // Add ids associated with the group
389 ids.insert(iter.val());
390 }
391 }
392 }
393
394 if (ids.empty())
395 {
396 return PtrListOps::findMatching(*this, matcher);
397 }
398 else
399 {
400 ids.insert(PtrListOps::findMatching(*this, matcher));
401 }
402 }
403 else
404 {
405 // Literal string.
406 // Special version of above for reduced memory footprint.
407
408 const label patchId = PtrListOps::firstMatching(*this, matcher);
409
410 if (patchId >= 0)
411 {
412 return labelList(one{}, patchId);
413 }
414 else if (checkGroups)
415 {
416 const auto iter = groupPatchIDs().cfind(matcher);
417
418 if (iter.good())
419 {
420 // Add ids associated with the group
421 ids.insert(iter.val());
422 }
424 }
425
426 return ids.sortedToc();
427}
428
429
431(
432 const wordRes& matcher,
433 const bool useGroups
434) const
435{
436 if (matcher.empty())
437 {
438 return labelList();
439 }
440 else if (matcher.size() == 1)
441 {
442 return this->indices(matcher.front(), useGroups);
443 }
444
445 labelHashSet ids;
446
447 // Only check groups if requested and they exist
448 if (useGroups && this->hasGroupIDs())
449 {
450 ids.reserve(this->size());
451
452 const auto& groupLookup = groupPatchIDs();
453 forAllConstIters(groupLookup, iter)
454 {
455 if (matcher(iter.key()))
456 {
457 // Add ids associated with the group
458 ids.insert(iter.val());
459 }
460 }
461 }
462
463 if (ids.empty())
464 {
465 return PtrListOps::findMatching(*this, matcher);
466 }
467 else
468 {
469 ids.insert(PtrListOps::findMatching(*this, matcher));
470 }
471
472 return ids.sortedToc();
473}
474
475
477(
478 const wordRes& allow,
479 const wordRes& deny,
480 const bool useGroups
481) const
482{
483 if (allow.empty() && deny.empty())
484 {
485 // Fast-path: select all
486 return identity(this->size());
487 }
488
489 const wordRes::filter matcher(allow, deny);
490
491 labelHashSet ids;
492
493 // Only check groups if requested and they exist
494 if (useGroups && this->hasGroupIDs())
495 {
496 ids.reserve(this->size());
497
498 const auto& groupLookup = groupPatchIDs();
499 forAllConstIters(groupLookup, iter)
500 {
501 if (matcher(iter.key()))
502 {
503 // Add ids associated with the group
504 ids.insert(iter.val());
505 }
506 }
507 }
508
509 if (ids.empty())
510 {
511 return PtrListOps::findMatching(*this, matcher);
512 }
513 else
514 {
515 ids.insert(PtrListOps::findMatching(*this, matcher));
516 }
517
518 return ids.sortedToc();
519}
520
521
523(
524 const word& patchName,
525 bool allowNotFound
526) const
527{
528 //return mesh().boundaryMesh().findPatchID(patchName);
529 if (patchName.empty())
530 {
531 return -1;
532 }
533
534 const label patchId = PtrListOps::firstMatching(*this, patchName);
535
536 if (patchId >= 0)
537 {
538 return patchId;
539 }
540
541 if (!allowNotFound)
542 {
544 << "Patch '" << patchName << "' not found. "
545 << "Available patch names";
546
547 if (polyMesh::defaultRegion != mesh_.name())
548 {
550 << " in region '" << mesh_.name() << "'";
551 }
552
554 << " include: " << names() << endl
555 << exit(FatalError);
556 }
557
558 // Patch not found
559 if (debug)
560 {
561 Pout<< "label pointBoundaryMesh::findPatchID(const word&) const"
562 << " Patch named " << patchName << " not found. "
563 << "Available patch names: " << names() << endl;
564 }
566 // Not found, return -1
567 return -1;
568}
569
570
571const Foam::pointPatch*
572Foam::pointBoundaryMesh::cfindPatch(const word& patchName) const
573{
574 const pointPatchList& patches = *this;
575
576 if (!patchName.empty())
577 {
578 // Note: get() handles out-of-range access properly
579 return patches.get(PtrListOps::firstMatching(patches, patchName));
580 }
581
582 return nullptr;
583}
584
585
586void Foam::pointBoundaryMesh::calcGeometry()
587{
589
590 if
591 (
592 pBufs.commsType() == Pstream::commsTypes::buffered
593 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
594 )
595 {
596 forAll(*this, patchi)
597 {
598 operator[](patchi).initGeometry(pBufs);
599 }
600
601 pBufs.finishedSends();
602
603 forAll(*this, patchi)
604 {
605 operator[](patchi).calcGeometry(pBufs);
606 }
607 }
608 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
609 {
610 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
611
612 // Dummy.
613 pBufs.finishedSends();
614
615 for (const auto& schedEval : patchSchedule)
616 {
617 const label patchi = schedEval.patch;
618
619 if (schedEval.init)
620 {
621 operator[](patchi).initGeometry(pBufs);
622 }
623 else
624 {
625 operator[](patchi).calcGeometry(pBufs);
626 }
627 }
628 }
629}
630
631
633{
635
636 if
637 (
638 pBufs.commsType() == Pstream::commsTypes::buffered
639 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
640 )
641 {
642 forAll(*this, patchi)
643 {
644 operator[](patchi).initMovePoints(pBufs, p);
645 }
646
647 pBufs.finishedSends();
648
649 forAll(*this, patchi)
650 {
651 operator[](patchi).movePoints(pBufs, p);
652 }
653 }
654 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
655 {
656 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
657
658 // Dummy.
659 pBufs.finishedSends();
660
661 for (const auto& schedEval : patchSchedule)
662 {
663 const label patchi = schedEval.patch;
664
665 if (schedEval.init)
666 {
667 operator[](patchi).initMovePoints(pBufs, p);
668 }
669 else
670 {
671 operator[](patchi).movePoints(pBufs, p);
672 }
673 }
674 }
675}
676
677
679{
681
682 if
683 (
684 pBufs.commsType() == Pstream::commsTypes::buffered
685 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
686 )
687 {
688 forAll(*this, patchi)
689 {
690 operator[](patchi).initUpdateMesh(pBufs);
691 }
692
693 pBufs.finishedSends();
694
695 forAll(*this, patchi)
696 {
697 operator[](patchi).updateMesh(pBufs);
698 }
699 }
700 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
701 {
702 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
703
704 // Dummy.
705 pBufs.finishedSends();
706
707 for (const auto& schedEval : patchSchedule)
708 {
709 const label patchi = schedEval.patch;
710
711 if (schedEval.init)
712 {
713 operator[](patchi).initUpdateMesh(pBufs);
714 }
715 else
716 {
717 operator[](patchi).updateMesh(pBufs);
718 }
719 }
720 }
721}
722
723
725(
726 const labelUList& oldToNew,
727 const bool validBoundary
728)
729{
730 // Change order of patches
731 pointPatchList::reorder(oldToNew);
732
733 // Adapt indices
734 pointPatchList& patches = *this;
735
736 forAll(patches, patchi)
737 {
738 patches[patchi].index() = patchi;
739 }
740
741 // Clear group-to-patch addressing. Note: could re-calculate
742 groupIDsPtr_.reset(nullptr);
743
744 if (validBoundary)
745 {
746 updateMesh();
747 }
748
749 if (debug)
750 {
751 pointPatchList& Patches = *this;
752
753 Pout<< "pointBoundaryMesh::reorder"
754 << "(const labelUList&, const bool): "
755 << "reordered pointBoundaryMesh:" << endl;
757 for (const auto& pp : Patches)
758 {
759 Pout<< indent
760 << "index:" << pp.index() << " patch:" << pp.name()
761 << " type:" << pp.type() << endl;
762 }
764 }
765}
766
767
769{
770 const pointPatchList& patches = *this;
771
773
774 forAll(patches, patchi)
775 {
776 os << indent << patches[patchi].name() << nl
778 << incrIndent << patches[patchi] << decrIndent
780 }
781
783
784 // Check state of IOstream
785 os.check("pointBoundaryMesh::writeData(Ostream& os) const");
786
787 return os.good();
788}
789
790
791//bool Foam::pointBoundaryMesh::writeObject
792//(
793// IOstreamOption
794//) const
795//{
796// return regIOobject::writeObject(fmt, ver, IOstream::UNCOMPRESSED);
797//}
798
799
800// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
bool empty() const noexcept
True if the hash table is empty.
Definition HashTable.H:353
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
bool isReadOptional() const noexcept
True if (LAZY_READ) bits are set [same as READ_IF_PRESENT].
bool isReadRequired() const noexcept
True if (MUST_READ | READ_MODIFIED) bits are set.
readOption readOpt() const noexcept
Get the read option.
writeOption writeOpt() const noexcept
Get the write option.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
fileName objectRelPath() const
The object path relative to the case.
Definition IOobject.C:581
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:45
bool good() const noexcept
True if next operation might succeed.
Definition IOstream.H:281
virtual const fileName & name() const
The name of the stream.
Definition IOstream.C:33
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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 the send phase as being finished.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void setSize(const label n)
Same as resize().
Definition PtrList.H:357
void resize_null(const label newLen)
Set the addressed list to the given size, deleting all existing entries. Afterwards the list contains...
Definition PtrListI.H:113
constexpr PtrList() noexcept
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
T & front()
Access first element of the list, position [0].
Definition UListI.H:239
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
@ scheduled
"scheduled" (MPI standard) : (MPI_Send, MPI_Recv)
Definition UPstream.H:83
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
static commsTypes defaultCommsType
Default commsType.
Definition UPstream.H:1045
const T & operator[](const label i) const
Return const reference to the element at given position. FatalError for bounds problem or nullptr....
Definition UPtrListI.H:289
void reorder(const labelUList &oldToNew, const bool check=false)
Reorder elements. Reordering must be unique (ie, shuffle).
Definition UPtrList.C:62
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
label count() const noexcept
The number of non-nullptr entries in the list.
Definition UPtrList.H:890
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< facePointPatch > New(const polyPatch &, const pointBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Type type(bool followLink=true, bool checkGzip=false) const
Return the directory entry type: UNDEFINED, FILE, DIRECTORY (or SYMLINK).
Definition fileName.C:353
const lduSchedule & patchSchedule() const noexcept
Order in which the patches should be initialised/evaluated corresponding to the schedule.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
A pointBoundaryMesh is a pointPatch list with registered IO, a reference to the associated pointMesh,...
label nNonProcessor() const
The number of patches before the first processor patch.
wordList physicalTypes() const
Return a list of physical types.
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
virtual bool writeData(Ostream &) const
writeData member function required by regIOobject
const pointMesh & mesh() const noexcept
Return the mesh reference.
void reorder(const labelUList &oldToNew, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
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.
labelList indices(const wordRe &matcher, const bool useGroups) const
The (sorted) patch indices for all matches, optionally matching patch groups.
void movePoints(const pointField &)
Correct pointBoundaryMesh after moving points.
friend class pointMesh
Declare friendship with pointMesh.
wordList names() const
Return a list of patch names.
void updateMesh()
Correct pointBoundaryMesh after topology update.
label nProcessorPatches() const
The number of processorPointPatch patches.
const pointPatch * cfindPatch(const word &patchName) const
Find patch by name and return const pointer.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
static autoPtr< pointPatch > New(const word &name, const dictionary &dict, const label index, const pointBoundaryMesh &)
Return a pointer to a new patch created on freestore. Returns null if not found.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
void close()
Close Istream.
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
bool headerOk()
Read and check header info. Does not check the headerClassName.
@ BEGIN_BLOCK
Begin block [isseparator].
Definition token.H:178
@ END_BLOCK
End block [isseparator].
Definition token.H:179
@ BEGIN_LIST
Begin list [isseparator].
Definition token.H:174
@ END_LIST
End list [isseparator].
Definition token.H:175
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
bool isPattern() const noexcept
The wordRe is a pattern, not a literal string.
Definition wordReI.H:104
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & names
label patchId(-1)
#define WarningInFunction
Report a warning using Foam::Warning.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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 handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition lduSchedule.H:46
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
PtrList< pointPatch > pointPatchList
Store lists of pointPatch as a PtrList.
Definition pointPatch.H:58
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Extract name (as a word) from an object, typically using its name() method.
Definition word.H:341
Extract type (as a word) from an object, typically using its type() method.
Definition word.H:362
Functor wrapper of allow/deny lists of wordRe for filtering.
Definition wordRes.H:275