48bool Foam::solidBodyFvGeometryScheme::markChanges
50 const pointField& oldPoints,
51 const pointField& currPoints,
52 bitSet& isChangedPoint,
53 bitSet& isChangedFace,
57 isChangedPoint.setSize(oldPoints.size());
60 forAll(isChangedPoint, pointi)
62 isChangedPoint.set(pointi, oldPoints[pointi] != currPoints[pointi]);
66 <<
"SBM --- Changed points:"
76 isChangedFace.setSize(
mesh_.nFaces());
77 isChangedFace =
false;
79 isChangedCell.setSize(
mesh_.nCells());
80 isChangedCell =
false;
82 const auto& pointFaces =
mesh_.pointFaces();
83 const auto& own =
mesh_.faceOwner();
84 const auto& nbr =
mesh_.faceNeighbour();
87 for (
const label pointi : isChangedPoint)
89 for (
const auto facei : pointFaces[pointi])
91 isChangedFace.set(facei);
93 isChangedCell.set(own[facei]);
94 if (facei <
mesh_.nInternalFaces())
96 isChangedCell.set(nbr[facei]);
102 <<
"SBM --- Changed cells:"
110void Foam::solidBodyFvGeometryScheme::setMeshMotionData()
112 if (!cacheInitialised_ || !cacheMotion_)
116 changedFaceIDs_.clear();
117 changedPatchIDs_.clear();
118 changedCellIDs_.clear();
120 const pointField& oldPoints = mesh_.oldPoints();
121 const pointField& currPoints = mesh_.points();
123 if (oldPoints.size() != currPoints.size())
126 <<
"Old and current points sizes must be the same. "
127 <<
"Old points:" << oldPoints.size()
128 <<
" Current points:" << currPoints.size()
187 const bool changed = markChanges
202 changedCellIDs_ = isChangedCell.toc();
209 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
211 if (isChangedFace[facei])
213 changedFaceIDs.append(facei);
214 changedPatchIDs.append(-1);
218 const auto&
pbm = mesh_.boundaryMesh();
219 for (label patchi = 0; patchi <
pbm.
size(); ++patchi)
223 for (
const label meshFacei :
pp.range())
225 if (isChangedFace[meshFacei])
227 changedFaceIDs.append(meshFacei);
228 changedPatchIDs.append(patchi);
233 changedFaceIDs_.
transfer(changedFaceIDs);
234 changedPatchIDs_.transfer(changedPatchIDs);
237 cacheInitialised_ =
true;
243Foam::solidBodyFvGeometryScheme::solidBodyFvGeometryScheme
250 partialUpdate_(
dict.getOrDefault<bool>(
"partialUpdate", true)),
251 cacheMotion_(
dict.getOrDefault<bool>(
"cacheMotion", true)),
252 cacheInitialised_(false),
258 <<
"partialUpdate:" << partialUpdate_
259 <<
" cacheMotion:" << cacheMotion_
272 mesh_.hasCellCentres()
273 && mesh_.hasFaceCentres()
274 && mesh_.hasCellVolumes()
275 && mesh_.hasFaceAreas();
280 <<
"Creating initial geometry using primitiveMesh::updateGeom"
293 const pointField& oldPoints = mesh_.oldPoints();
294 const pointField& currPoints = mesh_.points();
296 if (oldPoints.size() != currPoints.size())
299 <<
"Old and current points sizes must be the same. "
300 <<
"Old points:" << oldPoints.size()
301 <<
" Current points:" << currPoints.size()
306 const faceList& faces = mesh_.faces();
308 auto tmeshPhi(
const_cast<fvMesh&
>(mesh_).setPhi());
311 const scalar rdt = 1.0/mesh_.time().deltaTValue();
314 auto&
meshPhi = tmeshPhi.ref();
315 auto& meshPhii =
meshPhi.primitiveFieldRef();
316 auto& meshPhiBf =
meshPhi.boundaryFieldRef();
322 forAll(changedFaceIDs_, i)
324 const face&
f = faces[changedFaceIDs_[i]];
326 if (changedPatchIDs_[i] == -1)
328 const label facei = changedFaceIDs_[i];
329 meshPhii[facei] =
f.sweptVol(oldPoints, currPoints)*rdt;
333 const label patchi = changedPatchIDs_[i];
341 const label patchFacei = changedFaceIDs_[i] -
pp.start();
343 meshPhiBf[patchi][patchFacei] =
344 f.sweptVol(oldPoints, currPoints)*rdt;
349 if (partialUpdate_ && haveGeometry)
384 std::move(faceCentres),
385 std::move(faceAreas),
386 std::move(cellCentres),
387 std::move(cellVolumes)
392 for (
const auto&
p : mesh_.boundaryMesh())
395 <<
" sum(Sf)=" <<
sum(
p.faceAreas())
396 <<
" sum(mag(Sf))=" <<
sum(
mag(
p.faceAreas()))
404 <<
"Performing complete geometry clear and update" <<
endl;
426 cacheInitialised_ =
false;
442 (faceCentres.
size() != mesh_.nFaces())
443 || (faceAreas.
size() != mesh_.nFaces())
444 || (cellCentres.
size() != mesh_.nCells())
445 || (cellVolumes.
size() != mesh_.nCells())
465 bitSet isChangedPoint;
466 bitSet isChangedFace;
467 bitSet isChangedCell;
468 const bool changed = markChanges
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
const word & name() const noexcept
Return the object name.
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
void size(const label n)
Older name for setAddressableSize.
label size() const noexcept
The number of entries in the list.
Default geometry calculation scheme. Slight stabilisation for bad meshes.
virtual bool updateGeom(const pointField &points, const refPtr< pointField > &oldPoints, pointField &faceCentres, vectorField &faceAreas, pointField &cellCentres, scalarField &cellVolumes) const
Calculate geometry quantities using mesh topology and provided points. If oldPoints provided only doe...
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
labelList toc() const
The indices of the on bits as a sorted labelList.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A face is a list of labels corresponding to mesh vertices.
Abstract base class for geometry calculation schemes.
const fvMesh & mesh_
Hold reference to mesh.
const fvMesh & mesh() const
Return mesh reference.
Mesh data needed to do the Finite Volume discretisation.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
A patch is a list of labels that address the faces in the global face list.
void clearGeom()
Clear geometry.
virtual void updateGeom()
Update all geometric data.
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
A class for managing references or pointers (no reference counting).
Geometry calculation scheme that performs geometry updates only in regions where the mesh has changed...
virtual bool updateGeom(const pointField &points, const refPtr< pointField > &oldPoints, pointField &faceCentres, vectorField &faceAreas, pointField &cellCentres, scalarField &cellVolumes) const
Calculate geometry quantities using mesh topology and provided points. If oldPoints provided only doe...
virtual void movePoints()
Do what is necessary if the mesh has moved.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh for topology changes.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
List< face > faceList
List of faces.
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.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
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.
#define forAll(list, i)
Loop across all elements in list.