Loading...
Searching...
No Matches
GAMGAgglomeration.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-2016 OpenFOAM Foundation
9 Copyright (C) 2019-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 "GAMGAgglomeration.H"
30#include "lduMesh.H"
31#include "lduMatrix.H"
32#include "Time.H"
33#include "GAMGInterface.H"
36#include "IOmanip.H"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
46}
47
48
49// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50
52(
53 const label nCreatedLevels,
54 const bool doProcessorAgglomerate
55)
56{
57 nCells_.setSize(nCreatedLevels, 0);
58 restrictAddressing_.setSize(nCreatedLevels);
59 nFaces_.setSize(nCreatedLevels, 0);
60 faceRestrictAddressing_.setSize(nCreatedLevels);
61 faceFlipMap_.setSize(nCreatedLevels);
62 nPatchFaces_.setSize(nCreatedLevels);
63 patchFaceRestrictAddressing_.setSize(nCreatedLevels);
64 meshLevels_.setSize(nCreatedLevels);
65
66 // Have procCommunicator_ always, even if not procAgglomerating
67 procCommunicator_.setSize(nCreatedLevels + 1, -1);
68 if (doProcessorAgglomerate && processorAgglomerate())
69 {
70 procAgglomMap_.setSize(nCreatedLevels);
71 agglomProcIDs_.setSize(nCreatedLevels);
72 procAgglomCommunicator_.setSize(nCreatedLevels);
73 procCellOffsets_.setSize(nCreatedLevels);
74 procFaceMap_.setSize(nCreatedLevels);
75 procBoundaryMap_.setSize(nCreatedLevels);
76 procBoundaryFaceMap_.setSize(nCreatedLevels);
77
78 procAgglomeratorPtr_().agglomerate();
79 }
80}
81
82
84{
85 Info<< "GAMGAgglomeration:" << nl
86 << " local agglomerator : " << type() << nl;
87 if (processorAgglomerate())
88 {
89 Info<< " processor agglomerator : "
90 << procAgglomeratorPtr_().type() << nl
91 << nl;
92 }
93
94 Info<< setw(36) << "nCells"
95 << setw(20) << "nFaces/nCells"
96 << setw(20) << "nInterfaces"
97 << setw(20) << "nIntFaces/nCells"
98 << setw(12) << "profile"
99 << nl
100 << setw(8) << "Level"
101 << setw(8) << "nProcs"
102 << " "
103 << setw(8) << "avg"
104 << setw(8) << "max"
105 << " "
106 << setw(8) << "avg"
107 << setw(8) << "max"
108 << " "
109 << setw(8) << "avg"
110 << setw(8) << "max"
111 << " "
112 << setw(8) << "avg"
113 << setw(8) << "max"
114 //<< " "
115 << setw(12) << "avg"
116 << nl
117 << setw(8) << "-----"
118 << setw(8) << "------"
119 << " "
120 << setw(8) << "---"
121 << setw(8) << "---"
122 << " "
123 << setw(8) << "---"
124 << setw(8) << "---"
125 << " "
126 << setw(8) << "---"
127 << setw(8) << "---"
128 << " "
129 << setw(8) << "---"
130 << setw(8) << "---"
131 //<< " "
132 << setw(12) << "---"
133 //<< " "
134 << nl;
135
136 const label maxSize = returnReduce(size(), maxOp<label>());
137
138 for (label levelI = 0; levelI <= maxSize; levelI++)
139 {
140 label nProcs = 0;
141 label nCells = 0;
142 scalar faceCellRatio = 0;
143 label nInterfaces = 0;
144 label nIntFaces = 0;
145 scalar ratio = 0.0;
146 scalar profile = 0.0;
147
148 if (hasMeshLevel(levelI))
149 {
150 nProcs = 1;
151
152 const lduMesh& fineMesh = meshLevel(levelI);
153 nCells = fineMesh.lduAddr().size();
154 faceCellRatio =
155 scalar(fineMesh.lduAddr().lowerAddr().size())/nCells;
156
157 const lduInterfacePtrsList interfaces =
158 fineMesh.interfaces();
159 forAll(interfaces, i)
160 {
161 if (interfaces.set(i))
162 {
163 nInterfaces++;
164 nIntFaces += interfaces[i].faceCells().size();
165 }
166 }
167 ratio = scalar(nIntFaces)/nCells;
168
169 profile = fineMesh.lduAddr().band().second();
170 }
171
172 label totNprocs = returnReduce(nProcs, sumOp<label>());
173
174 label maxNCells = returnReduce(nCells, maxOp<label>());
175 label totNCells = returnReduce(nCells, sumOp<label>());
176
177 scalar maxFaceCellRatio =
178 returnReduce(faceCellRatio, maxOp<scalar>());
179 scalar totFaceCellRatio =
180 returnReduce(faceCellRatio, sumOp<scalar>());
181
182 label maxNInt = returnReduce(nInterfaces, maxOp<label>());
183 label totNInt = returnReduce(nInterfaces, sumOp<label>());
184
185 scalar maxRatio = returnReduce(ratio, maxOp<scalar>());
186 scalar totRatio = returnReduce(ratio, sumOp<scalar>());
187
188 scalar totProfile = returnReduce(profile, sumOp<scalar>());
189
190 const int oldPrecision = Info.stream().precision(4);
191
192 Info<< setw(8) << levelI
193 << setw(8) << totNprocs
194 << " "
195 << setw(8) << totNCells/totNprocs
196 << setw(8) << maxNCells
197 << " "
198 << setw(8) << totFaceCellRatio/totNprocs
199 << setw(8) << maxFaceCellRatio
200 << " "
201 << setw(8) << scalar(totNInt)/totNprocs
202 << setw(8) << maxNInt
203 << " "
204 << setw(8) << totRatio/totNprocs
205 << setw(8) << maxRatio
206 << setw(12) << totProfile/totNprocs
207 << nl;
209 Info.stream().precision(oldPrecision);
210 }
211 Info<< endl;
212}
213
214
216(
217 const label nCellsInCoarsestLevel,
218 const label nFineCells,
219 const label nCoarseCells,
220 const label comm
221) const
222{
223 const label nTotalCoarseCells =
224 returnReduce(nCoarseCells, sumOp<label>(), UPstream::msgType(), comm);
225 if (nTotalCoarseCells < Pstream::nProcs(comm)*nCellsInCoarsestLevel)
226 {
227 return false;
228 }
229 else
230 {
231 const label nTotalFineCells =
232 returnReduce(nFineCells, sumOp<label>(), UPstream::msgType(), comm);
233 return nTotalCoarseCells < nTotalFineCells;
234 }
235}
236
237
238// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
239
241(
242 const lduMesh& mesh,
244)
245:
246 MeshObject_type
247 (
248 controlDict.getOrDefault<word>
249 (
250 "name",
251 GAMGAgglomeration::typeName
252 ),
253 mesh
254 ),
255
256 maxLevels_(50),
257
258 updateInterval_(controlDict.getOrDefault<label>("updateInterval", 1)),
259 requireUpdate_(false),
260
261 nCellsInCoarsestLevel_
262 (
263 controlDict.getOrDefault<label>("nCellsInCoarsestLevel", 10)
264 ),
265 meshInterfaces_(mesh.interfaces()),
266 procAgglomeratorPtr_
267 (
268 (
269 (UPstream::nProcs(mesh.comm()) > 1)
270 && controlDict.found("processorAgglomerator")
271 )
272 ? GAMGProcAgglomeration::New
273 (
274 controlDict.get<word>("processorAgglomerator"),
275 *this,
277 )
278 : autoPtr<GAMGProcAgglomeration>()
279 ),
280
281 nCells_(maxLevels_),
282 restrictAddressing_(maxLevels_),
283 nFaces_(maxLevels_),
284 faceRestrictAddressing_(maxLevels_),
285 faceFlipMap_(maxLevels_),
286 nPatchFaces_(maxLevels_),
287 patchFaceRestrictAddressing_(maxLevels_),
288
289 meshLevels_(maxLevels_)
290{
291 // Limit the cells in the coarsest level based on the local number of
292 // cells. Note: 2 for pair-wise
294 max(1, min(mesh.lduAddr().size()/2, nCellsInCoarsestLevel_));
295
296 // Ensure all procs see the same nCellsInCoarsestLevel_
297 reduce(nCellsInCoarsestLevel_, minOp<label>());
298
301 {
302 procAgglomMap_.setSize(maxLevels_);
303 agglomProcIDs_.setSize(maxLevels_);
309 }
310}
311
312
314(
315 const lduMesh& mesh,
317)
318{
319
320
321
322 const GAMGAgglomeration* agglomPtr =
323 mesh.thisDb().cfindObject<GAMGAgglomeration>
324 (
325 controlDict.getOrDefault<word>
326 (
327 "name",
328 GAMGAgglomeration::typeName
329 )
330 );
331
332 if (agglomPtr)
333 {
334 if (agglomPtr->requireUpdate_)
335 {
336 mesh.thisDb().checkOut(const_cast<GAMGAgglomeration*>(agglomPtr));
338 }
339
340 return *agglomPtr;
341 }
342
343 {
344 const word agglomeratorType
345 (
346 controlDict.getOrDefault<word>("agglomerator", "faceAreaPair")
347 );
348
349 mesh.thisDb().time().libs().open
350 (
352 "geometricGAMGAgglomerationLibs",
353 lduMeshConstructorTablePtr_
354 );
355
356 auto* ctorPtr = lduMeshConstructorTable(agglomeratorType);
357
358 if (!ctorPtr)
359 {
361 << "Unknown GAMGAgglomeration type "
362 << agglomeratorType << ".\n"
363 << "Valid matrix GAMGAgglomeration types :"
364 << lduMatrixConstructorTablePtr_->sortedToc() << endl
365 << "Valid geometric GAMGAgglomeration types :"
366 << lduMeshConstructorTablePtr_->sortedToc()
367 << exit(FatalError);
368 }
369
370 auto agglomPtr(ctorPtr(mesh, controlDict));
371 if (debug)
372 {
373 agglomPtr().printLevels();
374 }
375 return regIOobject::store(agglomPtr);
376 }
377}
378
379
380const Foam::GAMGAgglomeration& Foam::GAMGAgglomeration::New
381(
382 const lduMatrix& matrix,
384)
385{
386 const lduMesh& mesh = matrix.mesh();
387
388 const GAMGAgglomeration* agglomPtr =
390 (
391 controlDict.getOrDefault<word>
392 (
393 "name",
394 GAMGAgglomeration::typeName
395 )
396 );
397
398 if (agglomPtr)
399 {
400 if (agglomPtr->requireUpdate_)
401 {
402 mesh.thisDb().checkOut(const_cast<GAMGAgglomeration*>(agglomPtr));
403 return GAMGAgglomeration::New(matrix, controlDict);
404 }
405
406 return *agglomPtr;
407 }
408
409 {
410 const word agglomeratorType
411 (
412 controlDict.getOrDefault<word>("agglomerator", "faceAreaPair")
413 );
414
415 mesh.thisDb().time().libs().open
416 (
418 "algebraicGAMGAgglomerationLibs",
419 lduMatrixConstructorTablePtr_
420 );
421
422 auto* ctorPtr = lduMatrixConstructorTable(agglomeratorType);
423
424 if (!ctorPtr)
425 {
426 return New(mesh, controlDict);
427 }
428 else
429 {
430 auto agglomPtr(ctorPtr(matrix, controlDict));
431 if (debug)
432 {
433 agglomPtr().printLevels();
435 return regIOobject::store(agglomPtr);
436 }
437 }
438}
439
440
441const Foam::GAMGAgglomeration& Foam::GAMGAgglomeration::New
442(
443 const lduMesh& mesh,
444 const scalarField& cellVolumes,
445 const vectorField& faceAreas,
447)
448{
449
450 const GAMGAgglomeration* agglomPtr =
452 (
453 controlDict.getOrDefault<word>
454 (
455 "name",
456 GAMGAgglomeration::typeName
457 )
458 );
459
460 if (agglomPtr)
461 {
462 if (agglomPtr->requireUpdate_)
463 {
464 mesh.thisDb().checkOut(const_cast<GAMGAgglomeration*>(agglomPtr));
466 }
467
468 return *agglomPtr;
469 }
470
471 {
472 const word agglomeratorType
473 (
474 controlDict.getOrDefault<word>("agglomerator", "faceAreaPair")
475 );
476
477 const_cast<Time&>(mesh.thisDb().time()).libs().open
478 (
480 "geometricGAMGAgglomerationLibs",
481 geometryConstructorTablePtr_
482 );
483
484 auto* ctorPtr = geometryConstructorTable(agglomeratorType);
485
486 if (!ctorPtr)
487 {
489 << "Unknown GAMGAgglomeration type "
490 << agglomeratorType << ".\n"
491 << "Valid geometric GAMGAgglomeration types :"
492 << geometryConstructorTablePtr_->sortedToc()
493 << exit(FatalError);
494 }
495
496 auto agglomPtr
497 (
498 ctorPtr
499 (
500 mesh,
501 cellVolumes,
502 faceAreas,
504 )
505 );
506 if (debug)
507 {
508 agglomPtr().printLevels();
509 }
510 return regIOobject::store(agglomPtr);
511 }
512}
513
514
515// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
518{}
519
520
521// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
522
525 return
526 (updateInterval_ > 0)
527 && ((mesh_.thisDb().time().timeIndex() % updateInterval_) == 0);
528}
529
530
532{
533 if (requiresUpdate())
534 {
536 }
537 // What to return?
538 return requireUpdate_;
539}
540
541
543(
544 const label i
545) const
546{
547 if (i == 0)
548 {
549 return mesh_;
550 }
551 else
552 {
553 return meshLevels_[i - 1];
554 }
555}
556
557
558bool Foam::GAMGAgglomeration::hasMeshLevel(const label i) const
559{
560 if (i == 0)
561 {
562 return true;
563 }
564 else
565 {
566 return meshLevels_.set(i - 1);
567 }
568}
569
570
572(
573 const label i
574) const
575{
576 if (i == 0)
577 {
578 return meshInterfaces_;
579 }
580 else
581 {
582 return meshLevels_[i - 1].rawInterfaces();
583 }
584}
585
586
587void Foam::GAMGAgglomeration::clearLevel(const label i)
588{
589 if (hasMeshLevel(i))
590 {
591 meshLevels_.set(i - 1, nullptr);
592
593 if (i < nCells_.size())
594 {
595 nCells_[i] = -555;
596 restrictAddressing_.set(i, nullptr);
597 nFaces_[i] = -666;
598 faceRestrictAddressing_.set(i, nullptr);
599 faceFlipMap_.set(i, nullptr);
600 nPatchFaces_.set(i, nullptr);
601 patchFaceRestrictAddressing_.set(i, nullptr);
602 }
603 }
604}
605
606
608(
609 const label leveli
610) const
611{
612 return procAgglomMap_[leveli];
613}
614
615
617(
618 const label leveli
619) const
620{
621 return agglomProcIDs_[leveli];
622}
623
625bool Foam::GAMGAgglomeration::hasProcMesh(const label leveli) const
626{
627 return procCommunicator_[leveli] != -1;
628}
629
631Foam::label Foam::GAMGAgglomeration::procCommunicator(const label leveli) const
632{
633 return procCommunicator_[leveli];
634}
635
637Foam::label Foam::GAMGAgglomeration::agglomCommunicator(const label leveli) const
638{
639 return procAgglomCommunicator_[leveli];
640}
641
642
644(
645 const label leveli
646) const
647{
648 return procCellOffsets_[leveli];
649}
650
651
653(
654 const label leveli
655) const
656{
657 return procFaceMap_[leveli];
658}
659
660
662(
663 const label leveli
664) const
665{
666 return procBoundaryMap_[leveli];
667}
668
669
671(
672 const label leveli
673) const
674{
675 return procBoundaryFaceMap_[leveli];
676}
677
678
680(
681 labelList& newRestrict,
682 label& nNewCoarse,
683 const lduAddressing& fineAddressing,
684 const labelUList& restriction,
685 const label nCoarse
686)
687{
688 if (fineAddressing.size() != restriction.size())
689 {
691 << "nCells:" << fineAddressing.size()
692 << " agglom:" << restriction.size()
693 << abort(FatalError);
694 }
695
696 // Seed (master) for every region
697 labelList master(identity(fineAddressing.size()));
698
699 // Now loop and transport master through region
700 const labelUList& lower = fineAddressing.lowerAddr();
701 const labelUList& upper = fineAddressing.upperAddr();
702
703 while (true)
704 {
705 label nChanged = 0;
706
707 forAll(lower, facei)
708 {
709 const label own = lower[facei];
710 const label nei = upper[facei];
711
712 if (restriction[own] == restriction[nei])
713 {
714 // coarse-mesh-internal face
715
716 if (master[own] < master[nei])
717 {
718 master[nei] = master[own];
719 nChanged++;
720 }
721 else if (master[own] > master[nei])
722 {
723 master[own] = master[nei];
724 nChanged++;
725 }
726 }
727 }
728
729 // reduce(nChanged, sumOp<label>());
730
731 if (nChanged == 0)
732 {
733 break;
734 }
735 }
736
737
738 // Count number of regions/masters per coarse cell
739 labelListList coarseToMasters(nCoarse);
740 nNewCoarse = 0;
741 forAll(restriction, celli)
742 {
743 labelList& masters = coarseToMasters[restriction[celli]];
744
745 if (!masters.found(master[celli]))
746 {
747 masters.append(master[celli]);
748 nNewCoarse++;
749 }
750 }
751
752 if (nNewCoarse > nCoarse)
753 {
754 //WarningInFunction
755 // << "Have " << nCoarse
756 // << " agglomerated cells but " << nNewCoarse
757 // << " disconnected regions" << endl;
758
759 // Keep coarseToMasters[0] the original coarse, allocate new ones
760 // for the others
761 labelListList coarseToNewCoarse(coarseToMasters.size());
762
763 nNewCoarse = nCoarse;
764
765 forAll(coarseToMasters, coarseI)
766 {
767 const labelList& masters = coarseToMasters[coarseI];
768
769 labelList& newCoarse = coarseToNewCoarse[coarseI];
770 newCoarse.setSize(masters.size());
771 newCoarse[0] = coarseI;
772 for (label i=1; i<newCoarse.size(); i++)
773 {
774 newCoarse[i] = nNewCoarse++;
775 }
776 }
777
778 newRestrict.setSize(fineAddressing.size());
779 forAll(restriction, celli)
780 {
781 const label coarseI = restriction[celli];
782
783 const label index = coarseToMasters[coarseI].find(master[celli]);
784 newRestrict[celli] = coarseToNewCoarse[coarseI][index];
785 }
786
787 return false;
788 }
789
790 return true;
791}
792
793
794// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
bool found
Geometric agglomerated algebraic multigrid agglomeration class.
const lduInterfacePtrsList meshInterfaces_
Cached mesh interfaces.
virtual bool movePoints()
Update when the mesh moves.
bool processorAgglomerate() const
Whether to agglomerate across processors.
void clearLevel(const label leveli)
const labelListListList & boundaryFaceMap(const label fineLeveli) const
Mapping from processor to procMesh boundary face.
label agglomCommunicator(const label fineLeveli) const
Communicator for collecting contributions.
static const GAMGAgglomeration & New(const lduMesh &mesh, const dictionary &controlDict)
Return the selected geometric agglomerator.
PtrList< labelListList > patchFaceRestrictAddressing_
Patch-local face restriction addressing array.
PtrList< boolList > faceFlipMap_
Face flip: for faces mapped to internal faces stores whether.
PtrList< UPstream::communicator > procAgglomCommunicator_
Communicator for collecting contributions. Note self-contained.
PtrList< labelList > nPatchFaces_
The number of (coarse) patch faces in each level.
const label maxLevels_
Max number of levels.
void compactLevels(const label nCreatedLevels, const bool doProcessorAgglomerate)
Shrink the number of levels to that specified. Optionally do.
static bool checkRestriction(labelList &newRestrict, label &nNewCoarse, const lduAddressing &fineAddressing, const labelUList &restriction, const label nCoarse)
Given restriction determines if coarse cells are connected.
label nCellsInCoarsestLevel_
Number of cells in coarsest level.
const labelList & procAgglomMap(const label fineLeveli) const
Mapping from processor to agglomerated processor (global, all processors have the same information)....
autoPtr< GAMGProcAgglomeration > procAgglomeratorPtr_
bool requireUpdate_
Does agglomeration require update.
labelList nFaces_
The number of (coarse) faces in each level.
PtrList< lduPrimitiveMesh > meshLevels_
Hierarchy of mesh addressing.
GAMGAgglomeration(const GAMGAgglomeration &)=delete
No copy construct.
const labelList & cellOffsets(const label fineLeveli) const
Mapping from processor to procMesh cells.
labelList nCells_
The number of cells in each level.
bool hasProcMesh(const label fineLeveli) const
Check that level has combined mesh.
bool continueAgglomerating(const label nCellsInCoarsestLevel, const label nCells, const label nCoarseCells, const label comm) const
Check the need for further agglomeration.
const labelListList & faceMap(const label fineLeveli) const
Mapping from processor to procMesh face.
bool requiresUpdate() const
Does the agglomeration need to be fully updated?
void printLevels() const
Print level overview.
labelList procCommunicator_
Communicator for given level.
PtrList< labelList > agglomProcIDs_
Per level the set of processors to agglomerate. Element 0 is.
PtrList< labelList > procAgglomMap_
Per level, per processor the processor it agglomerates into.
PtrList< labelListList > procBoundaryMap_
Mapping from processor to procMeshLevel boundary.
const labelListList & boundaryMap(const label fineLeveli) const
Mapping from processor to procMesh boundary.
const labelList & agglomProcIDs(const label fineLeveli) const
Set of processors to agglomerate. Element 0 is the master processor. (local, same only on those proce...
PtrList< labelList > procCellOffsets_
Mapping from processor to procMeshLevel cells.
PtrList< labelField > restrictAddressing_
Cell restriction addressing array.
PtrList< labelListList > procFaceMap_
Mapping from processor to procMeshLevel face.
PtrList< labelList > faceRestrictAddressing_
Face restriction addressing array.
friend class GAMGProcAgglomeration
Declare friendship with GAMGProcAgglomeration.
const lduInterfacePtrsList & interfaceLevel(const label leveli) const
Return LDU interface addressing of given level.
const label updateInterval_
Update agglomeration every updateInterval_ steps.
bool hasMeshLevel(const label leveli) const
Do we have mesh for given level?
label nCells(const label leveli) const
Return number of coarse cells (before processor agglomeration).
label procCommunicator(const label fineLeveli) const
Communicator for current level or -1.
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
PtrList< labelListListList > procBoundaryFaceMap_
Mapping from processor to procMeshLevel boundary face.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
dlLibraryTable & libs() const noexcept
Mutable access to the loaded dynamic libraries.
Definition Time.H:722
const T2 & second() const noexcept
Access the second element.
Definition Tuple2.H:142
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
Inter-processor communications stream.
Definition UPstream.H:69
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool open(bool verbose=true)
Open named, but unopened libraries. These names will normally have been added with push_back().
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Tuple2< label, scalar > band() const
Calculate bandwidth and profile of addressing.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
label size() const noexcept
Return number of equations.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition lduMatrix.H:753
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition lduMesh.H:54
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
virtual lduInterfacePtrsList interfaces() const =0
Return a list of pointers for each patch with only those pointing to interfaces being set.
const Time & time() const noexcept
Return time registry.
bool checkOut(regIOobject *io) const
Remove a regIOobject from registry and free memory if the object is ownedByRegistry....
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
bool store()
Register object with its registry and transfer ownership to the registry.
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
runTime controlDict().readEntry("adjustTimeStep"
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
List< labelListList > labelListListList
List of labelListList.
Definition labelList.H:41
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Omanip< int > setw(const int i)
Definition IOmanip.H:199
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
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...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
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
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299