dynamicOversetFvMesh.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) 2014-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 i
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
33#include "globalIndex.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
42}
43
44
45// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
46
48{
50
51 // The (processor-local part of the) stencil determines the local
52 // faces to add to the matrix. tbd: parallel
53 const labelListList& stencil = overlap.cellStencil();
54
55 // Get the base addressing
56 //const lduAddressing& baseAddr = dynamicMotionSolverFvMesh::lduAddr();
57 const lduAddressing& baseAddr = dynamicFvMesh::lduAddr();
58
59 // Add to the base addressing
60 labelList lowerAddr;
61 labelList upperAddr;
62 label nExtraFaces;
63
64 const globalIndex globalNumbering(baseAddr.size());
65 labelListList localFaceCells;
66 labelListList remoteFaceCells;
67
69 forAll(baseAddr, cellI)
70 {
71 globalCellIDs[cellI] = globalNumbering.toGlobal(cellI);
72 }
74
75
77 (
78 baseAddr,
79 stencil,
80 nExtraFaces,
81 lowerAddr,
82 upperAddr,
84 globalNumbering,
85 globalCellIDs,
86 localFaceCells,
87 remoteFaceCells
88 );
89
90 if (debug)
91 {
92 Pout<< "dynamicOversetFvMesh::update() : extended addressing from"
93 << " nFaces:" << baseAddr.lowerAddr().size()
94 << " to nFaces:" << lowerAddr.size()
95 << " nExtraFaces:" << nExtraFaces << endl;
96 }
97
98
99 // Now for the tricky bits. We want to hand out processor faces according
100 // to the localFaceCells/remoteFaceCells. Ultimately we need
101 // per entry in stencil:
102 // - the patch (or -1 for internal faces)
103 // - the face (is either an internal face index or a patch face index)
104
106
107 // Per processor to owner (local)/neighbour (remote)
109 List<DynamicList<label>> dynProcNeighbour(Pstream::nProcs());
110 forAll(stencil, celli)
111 {
112 const labelList& nbrs = stencil[celli];
113 stencilPatches_[celli].setSize(nbrs.size());
114 stencilPatches_[celli] = -1;
115
116 forAll(nbrs, nbri)
117 {
118 if (stencilFaces_[celli][nbri] == -1)
119 {
120 const label nbrCelli = nbrs[nbri];
121 label globalNbr = globalCellIDs[nbrCelli];
122 label proci = globalNumbering.whichProcID(globalNbr);
123 label remoteCelli = globalNumbering.toLocal(proci, globalNbr);
124
125 // Overwrite the face to be a patch face
126 stencilFaces_[celli][nbri] = procOwner[proci].size();
127 stencilPatches_[celli][nbri] = proci;
128 procOwner[proci].append(celli);
129 dynProcNeighbour[proci].append(remoteCelli);
130
131 //Pout<< "From neighbour proc:" << proci
132 // << " allocating patchFace:" << stencilFaces_[celli][nbri]
133 // << " to get remote cell " << remoteCelli
134 // << " onto local cell " << celli << endl;
135 }
136 }
137 }
138 labelListList procNeighbour(dynProcNeighbour.size());
139 forAll(procNeighbour, i)
140 {
141 procNeighbour[i] = std::move(dynProcNeighbour[i]);
142 }
143 labelListList mySendCells;
144 Pstream::exchange<labelList, label>(procNeighbour, mySendCells);
145
146 label nbri = 0;
147 forAll(procOwner, proci)
148 {
149 if (procOwner[proci].size())
150 {
151 nbri++;
152 }
153 if (mySendCells[proci].size())
154 {
155 nbri++;
156 }
157 }
158 remoteStencilInterfaces_.setSize(nbri);
159 nbri = 0;
160
161 // E.g. if proc1 needs some data from proc2 and proc2 needs some data from
162 // proc1. We first want the pair : proc1 receive and proc2 send
163 // and then the pair : proc1 send, proc2 receive
164
165
166 labelList procToInterface(Pstream::nProcs(), -1);
167
168 forAll(procOwner, proci)
169 {
170 if (proci < Pstream::myProcNo() && procOwner[proci].size())
171 {
172 if (debug)
173 {
174 Pout<< "Adding interface " << nbri
175 << " to receive my " << procOwner[proci].size()
176 << " from " << proci << endl;
177 }
178 procToInterface[proci] = nbri;
180 (
181 nbri++,
183 (
184 procOwner[proci],
186 proci,
187 tensorField(0),
189 )
190 );
191 }
192 else if (proci > Pstream::myProcNo() && mySendCells[proci].size())
193 {
194 if (debug)
195 {
196 Pout<< "Adding interface " << nbri
197 << " to send my " << mySendCells[proci].size()
198 << " to " << proci << endl;
199 }
201 (
202 nbri++,
204 (
205 mySendCells[proci],
207 proci,
208 tensorField(0),
210 )
211 );
212 }
213 }
214 forAll(procOwner, proci)
215 {
216 if (proci > Pstream::myProcNo() && procOwner[proci].size())
217 {
218 if (debug)
219 {
220 Pout<< "Adding interface " << nbri
221 << " to receive my " << procOwner[proci].size()
222 << " from " << proci << endl;
223 }
224 procToInterface[proci] = nbri;
226 (
227 nbri++,
229 (
230 procOwner[proci],
232 proci,
233 tensorField(0),
235 )
236 );
237 }
238 else if (proci < Pstream::myProcNo() && mySendCells[proci].size())
239 {
240 if (debug)
241 {
242 Pout<< "Adding interface " << nbri
243 << " to send my " << mySendCells[proci].size()
244 << " to " << proci << endl;
245 }
247 (
248 nbri++,
250 (
251 mySendCells[proci],
253 proci,
254 tensorField(0),
256 )
257 );
258 }
259 }
260
261
262 // Rewrite stencilPatches now we know the actual interface (procToInterface)
263 for (auto& patches : stencilPatches_)
264 {
265 for (auto& interface : patches)
266 {
267 if (interface != -1)
268 {
269 interface = procToInterface[interface]+boundary().size();
270 }
271 }
272 }
273
274
275 // Get addressing and interfaces of all interfaces
276
277
279 {
280 const fvBoundaryMesh& fvp = boundary();
281
282 patchAddr.setSize(fvp.size() + remoteStencilInterfaces_.size());
283
284 //allInterfaces_ = dynamicMotionSolverFvMesh::interfaces();
286 allInterfaces_.setSize(patchAddr.size());
287
288 forAll(fvp, patchi)
289 {
290 patchAddr.set(patchi, &fvp[patchi].faceCells());
291 }
293 {
294 const label patchi = fvp.size()+i;
297
298 //Pout<< "at patch:" << patchi
299 // << " have procPatch:" << pp.type()
300 // << " from:" << pp.myProcNo()
301 // << " to:" << pp.neighbProcNo()
302 // << " with fc:" << pp.faceCells().size() << endl;
303
304 patchAddr.set(patchi, &pp.faceCells());
305 allInterfaces_.set(patchi, &pp);
306 }
307 }
308 const lduSchedule ps
309 (
310 lduPrimitiveMesh::nonBlockingSchedule<processorLduInterface>
311 (
313 )
314 );
315
316 lduPtr_.reset
317 (
319 (
320 nCells(),
321 std::move(lowerAddr),
322 std::move(upperAddr),
323 patchAddr,
324 ps
325 )
326 );
327
328
329 // Check
330 if (debug)
331 {
332 const lduAddressing& addr = lduPtr_(); //this->lduAddr();
333
334 Pout<< "Adapted addressing:"
335 << " lower:" << addr.lowerAddr().size()
336 << " upper:" << addr.upperAddr().size() << endl;
337
338 // Using lduAddressing::patch
339 forAll(patchAddr, patchI)
340 {
341 Pout<< " " << patchI << "\tpatchAddr:"
342 << addr.patchAddr(patchI).size()
343 << endl;
344 }
345
346 // Using interfaces
347 const lduInterfacePtrsList& iFaces = allInterfaces_;
348 Pout<< "Adapted interFaces:" << iFaces.size() << endl;
349 forAll(iFaces, patchI)
350 {
351 if (iFaces.set(patchI))
352 {
353 Pout<< " " << patchI << "\tinterface:"
354 << iFaces[patchI].type() << endl;
355 }
356 }
357 }
358
359 return true;
360}
361
362
364(
365 const labelList& types,
366 const labelList& nbrTypes,
367 const scalarField& norm,
368 const scalarField& nbrNorm,
369 const label celli,
370 bitSet& isFront
371) const
372{
373 const labelList& own = faceOwner();
374 const labelList& nei = faceNeighbour();
375 const cell& cFaces = cells()[celli];
376
377 scalar avg = 0.0;
378 label n = 0;
379 label nFront = 0;
380 for (const label facei : cFaces)
381 {
382 if (isInternalFace(facei))
383 {
384 label nbrCelli = (own[facei] == celli ? nei[facei] : own[facei]);
385 if (norm[nbrCelli] == -GREAT)
386 {
387 // Invalid neighbour. Add to front
388 if (isFront.set(facei))
389 {
390 nFront++;
391 }
392 }
393 else
394 {
395 // Valid neighbour. Add to average
396 avg += norm[nbrCelli];
397 n++;
398 }
399 }
400 else
401 {
402 if (nbrNorm[facei-nInternalFaces()] == -GREAT)
403 {
404 if (isFront.set(facei))
405 {
406 nFront++;
407 }
408 }
409 else
410 {
411 avg += nbrNorm[facei-nInternalFaces()];
412 n++;
413 }
414 }
415 }
416
417 if (n > 0)
418 {
419 return avg/n;
420 }
421 else
422 {
423 return norm[celli];
424 }
425}
426
427
429(
430 const GAMGAgglomeration& agglom
431) const
432{
433 labelList cellToCoarse(identity(nCells()));
434 labelListList coarseToCell(invertOneToMany(nCells(), cellToCoarse));
435
436 // Write initial agglomeration
437 {
438 volScalarField scalarAgglomeration
439 (
441 (
442 "agglomeration",
443 this->time().timeName(),
444 *this,
447 false
448 ),
449 *this,
451 );
452 scalarField& fld = scalarAgglomeration.primitiveFieldRef();
453 forAll(fld, celli)
454 {
455 fld[celli] = cellToCoarse[celli];
456 }
457 fld /= max(fld);
459 <
462 >(scalarAgglomeration.boundaryFieldRef(), false);
463 scalarAgglomeration.write();
464
465 Info<< "Writing initial cell distribution to "
466 << this->time().timeName() << endl;
467 }
468
469
470 for (label level = 0; level < agglom.size(); level++)
471 {
472 const labelList& addr = agglom.restrictAddressing(level);
473 label coarseSize = max(addr)+1;
474
475 Info<< "Level : " << level << endl
476 << returnReduce(addr.size(), sumOp<label>()) << endl
477 << " current size : "
478 << returnReduce(addr.size(), sumOp<label>()) << endl
479 << " agglomerated size : "
480 << returnReduce(coarseSize, sumOp<label>()) << endl;
481
482 forAll(addr, fineI)
483 {
484 const labelList& cellLabels = coarseToCell[fineI];
485 forAll(cellLabels, i)
486 {
487 cellToCoarse[cellLabels[i]] = addr[fineI];
488 }
489 }
490 coarseToCell = invertOneToMany(coarseSize, cellToCoarse);
491
492 // Write agglomeration
493 {
494 volScalarField scalarAgglomeration
495 (
497 (
498 "agglomeration_" + Foam::name(level),
499 this->time().timeName(),
500 *this,
503 false
504 ),
505 *this,
507 );
508 scalarField& fld = scalarAgglomeration.primitiveFieldRef();
509 forAll(fld, celli)
510 {
511 fld[celli] = cellToCoarse[celli];
512 }
513 //if (normalise)
514 //{
515 // fld /= max(fld);
516 //}
518 <
521 >(scalarAgglomeration.boundaryFieldRef(), false);
522 scalarAgglomeration.write();
523 }
524 }
525}
526
527
528// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
529
531(
532 const IOobject& io,
533 const bool doInit
534)
535:
537{
538 if (doInit)
539 {
540 init(false); // do not initialise lower levels
541 }
542}
543
544
545bool Foam::dynamicOversetFvMesh::init(const bool doInit)
546{
547 if (doInit)
548 {
550 }
551
552 active_ = false;
553
554 // Load stencil (but do not update)
555 (void)Stencil::New(*this, false);
556
557 // Assume something changed
558 return true;
559}
560
561
562// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
563
565{}
566
567
568// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
569
571{
572 if (!active_)
573 {
574 //return dynamicMotionSolverFvMesh::lduAddr();
575 return dynamicFvMesh::lduAddr();
576 }
577 if (!lduPtr_)
578 {
579 // Build extended addressing
580 updateAddressing();
581 }
582 return *lduPtr_;
583}
584
585
587{
588 if (!active_)
589 {
590 //return dynamicMotionSolverFvMesh::interfaces();
592 }
593 if (!lduPtr_)
594 {
595 // Build extended addressing
596 updateAddressing();
597 }
598 return allInterfaces_;
599}
600
601
604{
605 if (!lduPtr_)
606 {
608 << "Extended addressing not allocated" << abort(FatalError);
609 }
610
611 return *lduPtr_;
612}
613
614
616{
617 //if (dynamicMotionSolverFvMesh::update())
619 {
620 // Calculate the local extra faces for the interpolation. Note: could
621 // let demand-driven lduAddr() trigger it but just to make sure.
622 updateAddressing();
623
624 // Addressing and/or weights have changed. Make interpolated cells
625 // up to date with donors
626 interpolateFields();
627
628 return true;
629 }
630
631 return false;
632}
633
634
636{
637 if (name.ends_with("_0"))
638 {
639 return baseName(name.substr(0, name.size()-2));
640 }
641
642 return name;
643}
644
645
647{
648 // Add the stencil suppression list
649 wordHashSet suppressed(Stencil::New(*this).nonInterpolatedFields());
650
651 // Use whatever the solver has set up as suppression list
652 const dictionary* dictPtr
653 (
654 this->schemesDict().findDict("oversetInterpolationSuppressed")
655 );
656 if (dictPtr)
657 {
658 suppressed.insert(dictPtr->toc());
659 }
660
661 interpolate<volScalarField>(suppressed);
662 interpolate<volVectorField>(suppressed);
663 interpolate<volSphericalTensorField>(suppressed);
664 interpolate<volSymmTensorField>(suppressed);
665 interpolate<volTensorField>(suppressed);
666
667 return true;
668}
669
670
671
673(
674 IOstreamOption streamOpt,
675 const bool valid
676) const
677{
678 //bool ok = dynamicMotionSolverFvMesh::writeObject(streamOpt, valid);
679 bool ok = dynamicMotionSolverListFvMesh::writeObject(streamOpt, valid);
680
681 // For postprocessing : write cellTypes and zoneID
682 {
684
686
687 volScalarField volTypes
688 (
690 (
691 "cellTypes",
692 this->time().timeName(),
693 *this,
696 false
697 ),
698 *this,
700 zeroGradientFvPatchScalarField::typeName
701 );
702
703 forAll(volTypes.internalField(), cellI)
704 {
705 volTypes[cellI] = cellTypes[cellI];
706 }
707 volTypes.correctBoundaryConditions();
708 volTypes.writeObject(streamOpt, valid);
709 }
710 {
711 volScalarField volZoneID
712 (
714 (
715 "zoneID",
716 this->time().timeName(),
717 *this,
720 false
721 ),
722 *this,
724 zeroGradientFvPatchScalarField::typeName
725 );
726
728 const labelIOList& zoneID = overlap.zoneID();
729
730 forAll(zoneID, cellI)
731 {
732 volZoneID[cellI] = zoneID[cellI];
733 }
734 volZoneID.correctBoundaryConditions();
735 volZoneID.writeObject(streamOpt, valid);
736 }
737 if (debug)
738 {
740 const labelIOList& zoneID = overlap.zoneID();
741 const labelListList& cellStencil = overlap.cellStencil();
742
743 // Get remote zones
744 labelList donorZoneID(zoneID);
746
747 // Get remote cellCentres
748 pointField cc(C());
750
751 volScalarField volDonorZoneID
752 (
754 (
755 "donorZoneID",
756 this->time().timeName(),
757 *this,
760 false
761 ),
762 *this,
763 dimensionedScalar("minOne", dimless, scalar(-1)),
764 zeroGradientFvPatchScalarField::typeName
765 );
766
767 forAll(cellStencil, cellI)
768 {
769 const labelList& stencil = cellStencil[cellI];
770 if (stencil.size())
771 {
772 volDonorZoneID[cellI] = donorZoneID[stencil[0]];
773 for (label i = 1; i < stencil.size(); i++)
774 {
775 if (donorZoneID[stencil[i]] != volDonorZoneID[cellI])
776 {
777 WarningInFunction << "Mixed donor meshes for cell "
778 << cellI << " at " << C()[cellI]
779 << " donors:" << UIndirectList<point>(cc, stencil)
780 << endl;
781 volDonorZoneID[cellI] = -2;
782 }
783 }
784 }
785 }
786 //- Do not correctBoundaryConditions since re-interpolates!
787 //volDonorZoneID.correctBoundaryConditions();
788 volDonorZoneID.writeObject(streamOpt, valid);
789 }
790
791 return ok;
792}
793
794
795// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Geometric agglomerated algebraic multigrid agglomeration class.
const labelField & restrictAddressing(const label leveli) const
Return cell restrict addressing of given level.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
The IOstreamOption is a simple container for options an IOstream can normally have.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
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
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
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
void setSize(const label n)
Alias for resize()
Definition: UPtrList.H:261
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:590
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
virtual const labelUList & cellTypes() const
Return the cell type list.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
static const labelIOList & zoneID(const fvMesh &)
Helper: get reference to registered zoneID. Loads volScalarField.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
static int debug
Debug switch.
Definition: data.H:77
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
wordList toc() const
Return the table of contents.
Definition: dictionary.C:602
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:81
Dynamic mesh able to handle multiple motion solvers. NOTE: If the word entry "solvers" is not found i...
virtual bool update()
Update the mesh for both mesh motion and topology change.
dynamicFvMesh with support for overset meshes.
void writeAgglomeration(const GAMGAgglomeration &agglom) const
Debug: dump agglomeration.
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
virtual bool interpolateFields()
Update fields when mesh is updated.
virtual ~dynamicOversetFvMesh()
Destructor.
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended)
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
PtrList< const lduPrimitiveProcessorInterface > remoteStencilInterfaces_
Added (processor)lduInterfaces for remote bits of stencil.
autoPtr< fvMeshPrimitiveLduAddressing > lduPtr_
Extended addressing (extended with local interpolation stencils)
labelListList stencilFaces_
Corresponding faces (in above lduPtr) to the stencil.
virtual bool updateAddressing() const
Calculate the extended lduAddressing.
labelListList stencilPatches_
Corresponding patches (in above lduPtr) to the stencil.
const fvMeshPrimitiveLduAddressing & primitiveLduAddr() const
Return extended ldu addressing.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
virtual bool update()
Update the mesh for both mesh motion and topology change.
lduInterfacePtrsList allInterfaces_
Interfaces for above mesh. Contains both original and.
scalar cellAverage(const labelList &types, const labelList &nbrTypes, const scalarField &norm, const scalarField &nbrNorm, const label celli, bitSet &isFront) const
Average norm of valid neighbours.
labelList reverseFaceMap_
From old to new face labels.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
lduInterfacePtrsList interfaces() const
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: faMesh.C:904
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
Foam::fvBoundaryMesh.
Variant of fvMeshLduAddressing that contains addressing instead of slices.
static labelList addAddressing(const lduAddressing &addr, const labelListList &nbrCells, label &nExtraFaces, labelList &lower, labelList &upper, labelListList &nbrCellFaces, const globalIndex &, const labelList &globalCellIDs, labelListList &localFaceCells, labelListList &remoteFaceCells)
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:260
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:325
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
label size() const
Return number of equations.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
Concrete implementation of processor interface. Used to temporarily store settings.
virtual const labelUList & faceCells() const
Return faceCell addressing.
label constructSize() const noexcept
Constructed data size.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
label nCells() const noexcept
Number of mesh cells.
int myProcNo() const noexcept
Return processor number.
virtual bool write(const bool valid=true) const
Write using setting from DB.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
bool ends_with(const std::string &s) const
True if string ends with the given suffix (cf. C++20)
Definition: string.H:309
const word & baseName() const
Return const access to the base name of the sub-model.
Definition: subModelBase.C:119
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
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const cellShapeList & cells
const cellCellStencilObject & overlap
Definition: correctPhi.H:57
const labelIOList & zoneID
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:114
cellMask correctBoundaryConditions()
const labelList & cellTypes
Definition: setCellMask.H:33
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333