Loading...
Searching...
No Matches
oversetFvMeshBase.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-2022,2024 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
28#include "oversetFvMeshBase.H"
33#include "globalIndex.H"
34#include "GAMGAgglomeration.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
40 defineTypeNameAndDebug(oversetFvMeshBase, 0);
41}
42
43
44// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
45
47{
49
50 // The (processor-local part of the) stencil determines the local
51 // faces to add to the matrix. tbd: parallel
52 const labelListList& stencil = overlap.cellStencil();
53
54 // Get the base addressing
55 //const lduAddressing& baseAddr = dynamicMotionSolverFvMesh::lduAddr();
56 const lduAddressing& baseAddr = mesh_.fvMesh::lduAddr();
57
58 // Add to the base addressing
59 labelList lowerAddr;
60 labelList upperAddr;
61 label nExtraFaces;
62
63 const globalIndex globalNumbering(baseAddr.size());
64 labelListList localFaceCells;
65 labelListList remoteFaceCells;
66
67 labelList globalCellIDs(overlap.cellInterpolationMap().constructSize());
68 forAll(baseAddr, cellI)
69 {
70 globalCellIDs[cellI] = globalNumbering.toGlobal(cellI);
71 }
72 overlap.cellInterpolationMap().distribute(globalCellIDs);
73
74
76 (
77 baseAddr,
78 stencil,
79 nExtraFaces,
80 lowerAddr,
81 upperAddr,
83 globalNumbering,
84 globalCellIDs,
85 localFaceCells,
86 remoteFaceCells
87 );
88
89 if (debug)
90 {
91 Pout<< "oversetFvMeshBase::update() : extended addressing from"
92 << " nFaces:" << baseAddr.lowerAddr().size()
93 << " to nFaces:" << lowerAddr.size()
94 << " nExtraFaces:" << nExtraFaces << endl;
95 }
96
97 // Now for the tricky bits. We want to hand out processor faces according
98 // to the localFaceCells/remoteFaceCells. Ultimately we need
99 // per entry in stencil:
100 // - the patch (or -1 for internal faces)
101 // - the face (is either an internal face index or a patch face index)
102
103 stencilPatches_.setSize(stencilFaces_.size());
104
105 // Per processor to owner (local)/neighbour (remote)
106 List<DynamicList<label>> procOwner(Pstream::nProcs());
107 List<DynamicList<label>> dynProcNeighbour(Pstream::nProcs());
108 forAll(stencil, celli)
109 {
110 const labelList& nbrs = stencil[celli];
111 stencilPatches_[celli].setSize(nbrs.size());
112 stencilPatches_[celli] = -1;
113
114 forAll(nbrs, nbri)
115 {
116 if (stencilFaces_[celli][nbri] == -1)
117 {
118 const label nbrCelli = nbrs[nbri];
119 label globalNbr = globalCellIDs[nbrCelli];
120 label proci = globalNumbering.whichProcID(globalNbr);
121 label remoteCelli = globalNumbering.toLocal(proci, globalNbr);
122
123 // Overwrite the face to be a patch face
124 stencilFaces_[celli][nbri] = procOwner[proci].size();
125 stencilPatches_[celli][nbri] = proci;
126 procOwner[proci].append(celli);
127 dynProcNeighbour[proci].append(remoteCelli);
128
129 //Pout<< "From neighbour proc:" << proci
130 // << " allocating patchFace:" << stencilFaces_[celli][nbri]
131 // << " to get remote cell " << remoteCelli
132 // << " onto local cell " << celli << endl;
133 }
134 }
135 }
136 labelListList procNeighbour(dynProcNeighbour.size());
137 forAll(procNeighbour, i)
138 {
139 procNeighbour[i] = std::move(dynProcNeighbour[i]);
140 }
141 labelListList mySendCells;
142 Pstream::exchange<labelList, label>(procNeighbour, mySendCells);
143
144 label nbri = 0;
145 forAll(procOwner, proci)
146 {
147 if (procOwner[proci].size())
148 {
149 nbri++;
150 }
151 if (mySendCells[proci].size())
152 {
153 nbri++;
154 }
155 }
156 remoteStencilInterfaces_.setSize(nbri);
157 nbri = 0;
158
159 // E.g. if proc1 needs some data from proc2 and proc2 needs some data from
160 // proc1. We first want the pair : proc1 receive and proc2 send
161 // and then the pair : proc1 send, proc2 receive
162
163 labelList procToInterface(Pstream::nProcs(), -1);
164
165 forAll(procOwner, proci)
166 {
167 if (proci < Pstream::myProcNo() && procOwner[proci].size())
168 {
169 if (debug)
170 {
171 Pout<< "Adding interface " << nbri
172 << " to receive my " << procOwner[proci].size()
173 << " from " << proci << endl;
174 }
175 procToInterface[proci] = nbri;
177 (
178 nbri++,
179 new lduPrimitiveProcessorInterface
180 (
181 procOwner[proci],
183 proci,
184 tensorField(0),
186 )
187 );
188 }
189 else if (proci > Pstream::myProcNo() && mySendCells[proci].size())
190 {
191 if (debug)
192 {
193 Pout<< "Adding interface " << nbri
194 << " to send my " << mySendCells[proci].size()
195 << " to " << proci << endl;
196 }
198 (
199 nbri++,
200 new lduPrimitiveProcessorInterface
201 (
202 mySendCells[proci],
204 proci,
205 tensorField(0),
207 )
208 );
209 }
210 }
211 forAll(procOwner, proci)
212 {
213 if (proci > Pstream::myProcNo() && procOwner[proci].size())
214 {
215 if (debug)
216 {
217 Pout<< "Adding interface " << nbri
218 << " to receive my " << procOwner[proci].size()
219 << " from " << proci << endl;
220 }
221 procToInterface[proci] = nbri;
223 (
224 nbri++,
225 new lduPrimitiveProcessorInterface
226 (
227 procOwner[proci],
229 proci,
230 tensorField(0),
232 )
233 );
234 }
235 else if (proci < Pstream::myProcNo() && mySendCells[proci].size())
236 {
237 if (debug)
238 {
239 Pout<< "Adding interface " << nbri
240 << " to send my " << mySendCells[proci].size()
241 << " to " << proci << endl;
242 }
244 (
245 nbri++,
246 new lduPrimitiveProcessorInterface
247 (
248 mySendCells[proci],
250 proci,
251 tensorField(0),
253 )
254 );
255 }
256 }
257
258
259 // Rewrite stencilPatches now we know the actual interface (procToInterface)
260 for (auto& patches : stencilPatches_)
261 {
262 for (auto& interface : patches)
263 {
264 if (interface != -1)
265 {
266 interface = procToInterface[interface]+mesh_.boundary().size();
267 }
268 }
269 }
270
271 // Get addressing and interfaces of all interfaces
272
273 UPtrList<const labelUList> patchAddr;
274 {
275 const fvBoundaryMesh& fvp = mesh_.boundary();
276
277 patchAddr.setSize(fvp.size() + remoteStencilInterfaces_.size());
278
279 //allInterfaces_ = dynamicMotionSolverFvMesh::interfaces();
280 allInterfaces_ = mesh_.fvMesh::interfaces();
281 allInterfaces_.setSize(patchAddr.size());
282
283 forAll(fvp, patchi)
284 {
285 patchAddr.set(patchi, &fvp[patchi].faceCells());
286 }
287
289 {
290 label patchi = fvp.size()+i;
291 const lduPrimitiveProcessorInterface& pp =
293
294 //Pout<< "at patch:" << patchi
295 // << " have procPatch:" << pp.type()
296 // << " from:" << pp.myProcNo()
297 // << " to:" << pp.neighbProcNo()
298 // << " with fc:" << pp.faceCells().size() << endl;
299
300 patchAddr.set(patchi, &pp.faceCells());
301 allInterfaces_.set(patchi, &pp);
302 }
303 }
304
305 const lduSchedule ps
306 (
308 (
310 )
311 );
312
313 lduPtr_.reset
314 (
315 new fvMeshPrimitiveLduAddressing
316 (
317 mesh_.nCells(),
318 std::move(lowerAddr),
319 std::move(upperAddr),
320 patchAddr,
321 ps
322 )
323 );
324
325
326 // Check
327 if (debug)
328 {
329 const lduAddressing& addr = lduPtr_(); //this->lduAddr();
330
331 Pout<< "Adapted addressing:"
332 << " lower:" << addr.lowerAddr().size()
333 << " upper:" << addr.upperAddr().size() << endl;
334
335 // Using lduAddressing::patch
336 forAll(patchAddr, patchi)
337 {
338 Pout<< " " << patchi << "\tpatchAddr:"
339 << addr.patchAddr(patchi).size()
340 << endl;
341 }
342
343 // Using interfaces
344 const lduInterfacePtrsList& iFaces = mesh_.interfaces();
345 Pout<< "Adapted interFaces:" << iFaces.size() << endl;
346 forAll(iFaces, patchi)
347 {
348 if (iFaces.set(patchi))
349 {
350 Pout<< " " << patchi << "\tinterface:"
351 << iFaces[patchi].type() << endl;
352 }
354 }
355
356 return true;
357}
358
359
361(
362 const labelUList& types,
363 const labelUList& nbrTypes,
364 const scalarField& norm,
365 const scalarField& nbrNorm,
366 const label celli,
367 bitSet& isFront
368) const
369{
370 const labelList& own = mesh_.faceOwner();
371 const labelList& nei = mesh_.faceNeighbour();
372 const cell& cFaces = mesh_.cells()[celli];
373
374 scalar avg = 0.0;
375 label nTotal = 0;
376 for (const label facei : cFaces)
377 {
378 if (mesh_.isInternalFace(facei))
379 {
380 label nbrCelli = (own[facei] == celli ? nei[facei] : own[facei]);
381 if (norm[nbrCelli] == -GREAT)
382 {
383 // Invalid neighbour. Add to front
384 isFront.set(facei);
385 }
386 else
387 {
388 // Valid neighbour. Add to average
389 avg += norm[nbrCelli];
390 ++nTotal;
391 }
392 }
393 else
394 {
395 if (nbrNorm[facei-mesh_.nInternalFaces()] == -GREAT)
396 {
397 isFront.set(facei);
398 }
399 else
400 {
401 avg += nbrNorm[facei-mesh_.nInternalFaces()];
402 ++nTotal;
403 }
404 }
405 }
406
407 if (nTotal)
408 {
409 return avg/nTotal;
410 }
411 else
412 {
413 return norm[celli];
414 }
415}
416
417
419(
420 const GAMGAgglomeration& agglom
421) const
422{
423 labelList cellToCoarse(identity(mesh_.nCells()));
424 labelListList coarseToCell(invertOneToMany(mesh_.nCells(), cellToCoarse));
425
426 // Write initial agglomeration
427 {
428 volScalarField scalarAgglomeration
429 (
430 IOobject
431 (
432 "agglomeration",
433 mesh_.time().timeName(),
434 mesh_,
438 ),
439 mesh_,
441 );
442 scalarField& fld = scalarAgglomeration.primitiveFieldRef();
443 forAll(fld, celli)
444 {
445 fld[celli] = cellToCoarse[celli];
446 }
447 fld /= max(fld);
449 <
452 false
453 >(scalarAgglomeration.boundaryFieldRef());
454 scalarAgglomeration.write();
455
456 Info<< "Writing initial cell distribution to "
457 << mesh_.time().timeName() << endl;
458 }
459
460
461 for (label level = 0; level < agglom.size(); level++)
462 {
463 const labelList& addr = agglom.restrictAddressing(level);
464 label coarseSize = max(addr)+1;
465
466 Info<< "Level : " << level << endl
467 << " current size : "
468 << returnReduce(addr.size(), sumOp<label>()) << endl
469 << " agglomerated size : "
470 << returnReduce(coarseSize, sumOp<label>()) << endl;
471
472 forAll(addr, fineI)
473 {
474 const labelList& cellLabels = coarseToCell[fineI];
475 forAll(cellLabels, i)
476 {
477 cellToCoarse[cellLabels[i]] = addr[fineI];
478 }
479 }
480 coarseToCell = invertOneToMany(coarseSize, cellToCoarse);
481
482 // Write agglomeration
483 {
484 volScalarField scalarAgglomeration
485 (
487 (
488 "agglomeration_" + Foam::name(level),
489 mesh_.time().timeName(),
490 mesh_,
494 ),
495 mesh_,
497 );
498 scalarField& fld = scalarAgglomeration.primitiveFieldRef();
499 forAll(fld, celli)
500 {
501 fld[celli] = cellToCoarse[celli];
502 }
503 //if (normalise)
504 //{
505 // fld /= max(fld);
506 //}
508 <
511 false
512 >(scalarAgglomeration.boundaryFieldRef());
513 scalarAgglomeration.write();
515 }
516}
517
518
519// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
520
521Foam::oversetFvMeshBase::oversetFvMeshBase(const fvMesh& mesh, bool doInit)
522:
523 mesh_(mesh),
524 active_(false)
525{
526 // Load stencil (but do not update)
527 (void)Stencil::New(mesh_, false);
528}
529
530
531// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
534{}
535
536
537// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
538
540{
541 if (!active_)
542 {
543 //return dynamicMotionSolverFvMesh::lduAddr();
544 return mesh_.fvMesh::lduAddr();
545 }
546 if (!lduPtr_)
547 {
548 // Build extended addressing
550 }
551 return *lduPtr_;
552}
553
554
556{
557 if (!active_)
558 {
559 return mesh_.fvMesh::interfaces();
560 }
561 if (!lduPtr_)
562 {
563 // Build extended addressing
565 }
566 return allInterfaces_;
567}
568
569
572{
573 if (!lduPtr_)
574 {
576 << "Extended addressing not allocated" << abort(FatalError);
577 }
578
579 return *lduPtr_;
580}
581
582
585 // Cell-cell stencil is already mesh object. Clear out local
586 // addressing to force rebuilding addressing
587 lduPtr_.clear();
588}
589
590
592{
593 {
594 // Calculate the local extra faces for the interpolation. Note: could
595 // let demand-driven lduAddr() trigger it but just to make sure.
596 updateAddressing();
597
598 // Addressing and/or weights have changed. Make interpolated cells
599 // up to date with donors
600 interpolateFields();
601
602 return true;
603 }
604
605 return false;
606}
607
608
610{
611 const cellCellStencilObject& overlap = Stencil::New(mesh_);
612
613 // Add the stencil suppression list
614 wordHashSet suppressed(overlap.nonInterpolatedFields());
615
616 // Use whatever the solver has set up as suppression list
617 const dictionary* dictPtr
618 (
619 mesh_.schemesDict().findDict("oversetInterpolationSuppressed")
620 );
621 if (dictPtr)
622 {
623 suppressed.insert(dictPtr->toc());
624 }
625
626 overlap.interpolate<volScalarField>(mesh_, suppressed);
627 overlap.interpolate<volVectorField>(mesh_, suppressed);
628 overlap.interpolate<volSphericalTensorField>(mesh_, suppressed);
629 overlap.interpolate<volSymmTensorField>(mesh_, suppressed);
630 overlap.interpolate<volTensorField>(mesh_, suppressed);
631
632 return true;
633}
634
635
637(
638 IOstreamOption streamOpt,
639 const bool writeOnProc
640) const
641{
642 // For postprocessing : write cellTypes and zoneID
643 bool ok = true;
644 {
646
647 const labelUList& cellTypes = overlap.cellTypes();
648
649 volScalarField volTypes
650 (
652 (
653 "cellTypes",
654 mesh_.time().timeName(),
655 mesh_,
659 ),
660 mesh_,
663 );
664
665 forAll(volTypes.internalField(), cellI)
666 {
667 volTypes[cellI] = cellTypes[cellI];
668 }
669 volTypes.correctBoundaryConditions();
670 volTypes.writeObject(streamOpt, writeOnProc);
671 }
672 {
673 volScalarField volZoneID
674 (
675 IOobject
676 (
677 "zoneID",
678 mesh_.time().timeName(),
679 mesh_,
683 ),
684 mesh_,
687 );
688
689 const cellCellStencilObject& overlap = Stencil::New(mesh_);
690 const labelIOList& zoneID = overlap.zoneID();
691
692 forAll(zoneID, cellI)
693 {
694 volZoneID[cellI] = zoneID[cellI];
695 }
696 volZoneID.correctBoundaryConditions();
697 volZoneID.writeObject(streamOpt, writeOnProc);
698 }
699 if (debug)
700 {
702 const labelIOList& zoneID = overlap.zoneID();
703 const labelListList& cellStencil = overlap.cellStencil();
704
705 // Get remote zones
706 labelList donorZoneID(zoneID);
708
709 // Get remote cellCentres
710 pointField cc(mesh_.C());
712
713 volScalarField volDonorZoneID
714 (
716 (
717 "donorZoneID",
718 mesh_.time().timeName(),
719 mesh_,
723 ),
724 mesh_,
725 scalar(-1),
726 dimless,
728 );
729
730 forAll(cellStencil, cellI)
731 {
732 const labelList& stencil = cellStencil[cellI];
733 if (stencil.size())
734 {
735 volDonorZoneID[cellI] = donorZoneID[stencil[0]];
736 for (label i = 1; i < stencil.size(); i++)
737 {
738 if (donorZoneID[stencil[i]] != volDonorZoneID[cellI])
739 {
740 WarningInFunction << "Mixed donor meshes for cell "
741 << cellI << " at " << mesh_.C()[cellI]
742 << " donors:" << UIndirectList<point>(cc, stencil)
743 << endl;
744 volDonorZoneID[cellI] = -2;
745 }
746 }
747 }
748 }
749 //- Do not correctBoundaryConditions since re-interpolates!
750 //volDonorZoneID.correctBoundaryConditions();
752 <
755 >
756 (
757 volDonorZoneID
758 );
759 ok = volDonorZoneID.writeObject(streamOpt, writeOnProc);
760 }
761
762 return ok;
763}
764
765
766// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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 values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
void correctBoundaryConditions()
Correct boundary field.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
@ NO_REGISTER
Do not request registration (bool: false).
@ 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
A simple container for options an IOstream can normally have.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
static FOAM_NO_DANGLING_REFERENCE const cellCellStencilObject & New(const fvMesh &mesh, Args &&... args)
static void exchange(const UList< Container > &sendBufs, const labelUList &recvSizes, List< Container > &recvBufs, const int tag=UPstream::msgType(), const int comm=UPstream::worldComm, const bool wait=true)
Helper: exchange contiguous data. Sends sendBufs, receives into recvBufs using predetermined receive ...
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
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
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
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
void setSize(const label n)
Alias for resize().
Definition UPtrList.H:842
label size() const noexcept
The number of entries 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:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
static void correctBoundaryConditions(GeoField &psi)
Version of correctBoundaryConditions that excludes 'overset' bcs.
static void interpolate(Field< T > &psi, const fvMesh &mesh, const cellCellStencil &overlap, const List< scalarList > &wghts)
Interpolation of acceptor cells from donor cells.
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:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
wordList toc() const
Return the table of contents.
Definition dictionary.C:587
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
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)
Given additional addressing (in the form of additional neighbour cells, i.e. like cellCells).
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label toLocal(const label proci, const label i) const
From global to local on proci.
label whichProcID(const label proci, const label i) const
Which processor does global id come from? Checks proci first (assumed to occur reasonably frequently)...
label toGlobal(const label proci, const label i) const
From local to global on proci.
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
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.
label size() const noexcept
Return number of equations.
static lduSchedule nonBlockingSchedule(const lduInterfacePtrsList &)
Get non-scheduled send/receive schedule.
Concrete implementation of processor interface. Used to temporarily store settings.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
Support for overset functionality.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
bool active_
Select base addressing (false) or locally stored extended.
void writeAgglomeration(const GAMGAgglomeration &agglom) const
Debug: dump agglomeration.
static void correctBoundaryConditions(typename GeoField::Boundary &bfld)
Correct boundary conditions of certain type (TypeOnly = true) or explicitly not of the type (TypeOnly...
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
virtual bool interpolateFields()
Update fields when mesh is updated.
const fvMesh & mesh_
Reference to mesh.
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended).
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 update()
Update the mesh for both mesh motion and topology change.
lduInterfacePtrsList allInterfaces_
Interfaces for above mesh. Contains both original and above added processorLduInterfaces.
virtual ~oversetFvMeshBase()
Destructor.
scalar cellAverage(const labelUList &types, const labelUList &nbrTypes, const scalarField &norm, const scalarField &nbrNorm, const label celli, bitSet &isFront) const
Average norm of valid neighbours.
void clearOut()
Clear out local storage.
labelList reverseFaceMap_
From old to new face labels.
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
virtual void write(Ostream &os) const
Write.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const cellCellStencilObject & overlap
Definition correctPhi.H:57
#define WarningInFunction
Report a warning using Foam::Warning.
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
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition lduSchedule.H:46
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
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
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition ListOps.C:125
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const labelUList & cellTypes
Definition setCellMask.H:27
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299