Loading...
Searching...
No Matches
solidBodyFvGeometryScheme.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) 2021-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
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
30#include "surfaceFields.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
43 dict
44 );
45}
46
47
48bool Foam::solidBodyFvGeometryScheme::markChanges
49(
50 const pointField& oldPoints,
51 const pointField& currPoints,
52 bitSet& isChangedPoint,
53 bitSet& isChangedFace,
54 bitSet& isChangedCell
55) const
56{
57 isChangedPoint.setSize(oldPoints.size());
58
59 // Check for non-identical points
60 forAll(isChangedPoint, pointi)
61 {
62 isChangedPoint.set(pointi, oldPoints[pointi] != currPoints[pointi]);
63 }
64
66 << "SBM --- Changed points:"
67 << returnReduce(isChangedPoint.count(), sumOp<label>())
68 << endl;
69
70 // Quick return if no points have moved
71 if (returnReduceAnd(isChangedPoint.none()))
72 {
73 return false;
74 }
75
76 isChangedFace.setSize(mesh_.nFaces());
77 isChangedFace = false;
78
79 isChangedCell.setSize(mesh_.nCells());
80 isChangedCell = false;
81
82 const auto& pointFaces = mesh_.pointFaces();
83 const auto& own = mesh_.faceOwner();
84 const auto& nbr = mesh_.faceNeighbour();
85
86 // Identify faces and cells attached to moving points
87 for (const label pointi : isChangedPoint)
88 {
89 for (const auto facei : pointFaces[pointi])
90 {
91 isChangedFace.set(facei);
92
93 isChangedCell.set(own[facei]);
94 if (facei < mesh_.nInternalFaces())
95 {
96 isChangedCell.set(nbr[facei]);
97 }
98 }
99 }
100
102 << "SBM --- Changed cells:"
103 << returnReduce(isChangedCell.count(), sumOp<label>())
104 << endl;
105
106 return true;
107}
108
109
110void Foam::solidBodyFvGeometryScheme::setMeshMotionData()
111{
112 if (!cacheInitialised_ || !cacheMotion_)
113 {
114 DebugInFunction << "Creating cache" << endl;
115
116 changedFaceIDs_.clear(); // used for face areas, meshPhi
117 changedPatchIDs_.clear(); // used for meshPhi
118 changedCellIDs_.clear(); // used for cell volumes
119
120 const pointField& oldPoints = mesh_.oldPoints();
121 const pointField& currPoints = mesh_.points();
122
123 if (oldPoints.size() != currPoints.size())
124 {
126 << "Old and current points sizes must be the same. "
127 << "Old points:" << oldPoints.size()
128 << " Current points:" << currPoints.size()
129 << abort(FatalError);
130 }
131
132 //bitSet changedPoints(oldPoints.size());
133 //
135 //forAll(changedPoints, pointi)
136 //{
137 // changedPoints.set
138 // (
139 // pointi,
140 // oldPoints[pointi] != currPoints[pointi]
141 // );
142 //}
143 //
144 //DebugInfo
145 // << "SBM --- Changed points:"
146 // << returnReduce(changedPoints.count(), sumOp<label>())
147 // << endl;
148 //
150 //if (returnReduceAnd(changedPoints.none()))
151 //{
152 // return;
153 //}
154 //
155 //bitSet cellIDs(mesh_.nCells());
156 //bitSet faceIDs(mesh_.nFaces());
157 //
158 //const auto& pointFaces = mesh_.pointFaces();
159 //const auto& own = mesh_.faceOwner();
160 //const auto& nbr = mesh_.faceNeighbour();
161 //
163 //for (const label pointi : changedPoints)
164 //{
165 // for (const auto facei : pointFaces[pointi])
166 // {
167 // faceIDs.set(facei);
168 //
169 // cellIDs.set(own[facei]);
170 // if (facei < mesh_.nInternalFaces())
171 // {
172 // cellIDs.set(nbr[facei]);
173 // }
174 // }
175 //}
176 //
177 //changedCellIDs_ = cellIDs.toc();
178 //
179 //DebugInfo
180 // << "SBM --- Changed cells:"
181 // << returnReduce(changedCellIDs_.size(), sumOp<label>())
182 // << endl;
183
184 bitSet isChangedPoint;
185 bitSet isChangedFace;
186 bitSet isChangedCell;
187 const bool changed = markChanges
188 (
189 oldPoints,
190 currPoints,
191 isChangedPoint,
192 isChangedFace,
193 isChangedCell
194 );
195
196 // Quick return if no points have moved
197 if (!changed)
198 {
199 return;
200 }
201
202 changedCellIDs_ = isChangedCell.toc();
203
204
205 // Construct face and patch ID info
206
207 DynamicList<label> changedFaceIDs(isChangedFace.count());
208 DynamicList<label> changedPatchIDs(changedFaceIDs.capacity());
209 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
210 {
211 if (isChangedFace[facei])
212 {
213 changedFaceIDs.append(facei);
214 changedPatchIDs.append(-1);
215 }
216 }
217
218 const auto& pbm = mesh_.boundaryMesh();
219 for (label patchi = 0; patchi < pbm.size(); ++patchi)
220 {
221 const polyPatch& pp = pbm[patchi];
222
223 for (const label meshFacei : pp.range())
224 {
225 if (isChangedFace[meshFacei])
226 {
227 changedFaceIDs.append(meshFacei);
228 changedPatchIDs.append(patchi);
229 }
230 }
231 }
232
233 changedFaceIDs_.transfer(changedFaceIDs);
234 changedPatchIDs_.transfer(changedPatchIDs);
235 }
237 cacheInitialised_ = true;
238}
239
240
241// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
242
243Foam::solidBodyFvGeometryScheme::solidBodyFvGeometryScheme
244(
245 const fvMesh& mesh,
246 const dictionary& dict
247)
248:
250 partialUpdate_(dict.getOrDefault<bool>("partialUpdate", true)),
251 cacheMotion_(dict.getOrDefault<bool>("cacheMotion", true)),
252 cacheInitialised_(false),
253 changedFaceIDs_(),
254 changedPatchIDs_(),
255 changedCellIDs_()
256{
258 << "partialUpdate:" << partialUpdate_
259 << " cacheMotion:" << cacheMotion_
260 << endl;
261}
262
263
264// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
265
267{
268 // Note: not calling fvGeometryScheme::movePoints since we want to perform
269 // our own geometry manipulations
270
271 bool haveGeometry =
272 mesh_.hasCellCentres()
273 && mesh_.hasFaceCentres()
274 && mesh_.hasCellVolumes()
275 && mesh_.hasFaceAreas();
276
277 if (!haveGeometry)
278 {
280 << "Creating initial geometry using primitiveMesh::updateGeom"
281 << endl;
282
283 const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
284 return;
285 }
286
287 if (mesh_.moving())
288 {
289 setMeshMotionData();
290
291 DebugInFunction << "Performing partial meshPhi construction" << endl;
292
293 const pointField& oldPoints = mesh_.oldPoints();
294 const pointField& currPoints = mesh_.points();
295
296 if (oldPoints.size() != currPoints.size())
297 {
299 << "Old and current points sizes must be the same. "
300 << "Old points:" << oldPoints.size()
301 << " Current points:" << currPoints.size()
302 << abort(FatalError);
303 }
304
305 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
306 const faceList& faces = mesh_.faces();
307
308 auto tmeshPhi(const_cast<fvMesh&>(mesh_).setPhi());
309 if (tmeshPhi)
310 {
311 const scalar rdt = 1.0/mesh_.time().deltaTValue();
312
313 // Set the mesh flux
314 auto& meshPhi = tmeshPhi.ref();
315 auto& meshPhii = meshPhi.primitiveFieldRef();
316 auto& meshPhiBf = meshPhi.boundaryFieldRef();
317
318 //meshPhi == dimensionedScalar(dimVolume/dimTime, Zero);
319 meshPhii = Zero;
320 meshPhiBf == Zero;
321
322 forAll(changedFaceIDs_, i)
323 {
324 const face& f = faces[changedFaceIDs_[i]];
325
326 if (changedPatchIDs_[i] == -1)
327 {
328 const label facei = changedFaceIDs_[i];
329 meshPhii[facei] = f.sweptVol(oldPoints, currPoints)*rdt;
330 }
331 else
332 {
333 const label patchi = changedPatchIDs_[i];
334 const polyPatch& pp = pbm[patchi];
335
337 {
338 continue;
339 }
340
341 const label patchFacei = changedFaceIDs_[i] - pp.start();
342
343 meshPhiBf[patchi][patchFacei] =
344 f.sweptVol(oldPoints, currPoints)*rdt;
345 }
346 }
347 }
348
349 if (partialUpdate_ && haveGeometry)
350 {
351 // Keep base geometry and update as needed
352 DebugInFunction << "Performing partial geometry update" << endl;
353
354 // Initialise geometry using the old/existing values
355 vectorField faceCentres(mesh_.faceCentres());
356 vectorField faceAreas(mesh_.faceAreas());
357
358 // Make face centres and areas consistent with new points
360 (
361 mesh_,
362 changedFaceIDs_,
363 mesh_.points(),
364 faceCentres,
365 faceAreas
366 );
367
368 vectorField cellCentres(mesh_.cellCentres());
369 scalarField cellVolumes(mesh_.cellVolumes());
370
372 (
373 mesh_,
374 faceCentres,
375 faceAreas,
376 changedCellIDs_,
377 mesh_.cells(),
378 cellCentres,
379 cellVolumes
380 );
381
382 const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
383 (
384 std::move(faceCentres),
385 std::move(faceAreas),
386 std::move(cellCentres),
387 std::move(cellVolumes)
388 );
389
390 if (debug)
391 {
392 for (const auto& p : mesh_.boundaryMesh())
393 {
394 Pout<< "SBM --- " << p.name()
395 << " sum(Sf)=" << sum(p.faceAreas())
396 << " sum(mag(Sf))=" << sum(mag(p.faceAreas()))
397 << endl;
398 }
399 }
400 }
401 else
402 {
404 << "Performing complete geometry clear and update" << endl;
405
406 // Clear out old geometry
407 // Note: this recreates the old primitiveMesh::movePoints behaviour
408 const_cast<fvMesh&>(mesh_).primitiveMesh::clearGeom();
409
410 // Use lower level to calculate the geometry
411 const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
412 }
413 }
414 else
415 {
416 DebugInFunction << "Performing complete geometry update" << endl;
418 // Use lower level to calculate the geometry
419 const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
420 }
421}
422
425{
426 cacheInitialised_ = false;
427}
428
429
431(
432 const pointField& points,
433 const refPtr<pointField>& oldPoints,
434 pointField& faceCentres,
435 vectorField& faceAreas,
436 pointField& cellCentres,
437 scalarField& cellVolumes
438) const
439{
440 if
441 (
442 (faceCentres.size() != mesh_.nFaces())
443 || (faceAreas.size() != mesh_.nFaces())
444 || (cellCentres.size() != mesh_.nCells())
445 || (cellVolumes.size() != mesh_.nCells())
446 || !oldPoints
447 )
448 {
449 // Do all
451 (
452 points,
453 oldPoints,
454 faceCentres,
455 faceAreas,
456 cellCentres,
457 cellVolumes
458 );
459 }
460 else
461 {
462 // Since oldPoints provided assume that face & cell geometry is
463 // up to date with it
464
465 bitSet isChangedPoint;
466 bitSet isChangedFace;
467 bitSet isChangedCell;
468 const bool changed = markChanges
469 (
470 oldPoints(),
471 points,
472 isChangedPoint,
473 isChangedFace,
474 isChangedCell
475 );
476
477 if (!changed)
478 {
479 return false;
480 }
481
482 // Make face centres and areas consistent with new points
484 (
485 mesh_,
486 isChangedFace.toc(),
487 points,
488 faceCentres,
489 faceAreas
490 );
491
493 (
494 mesh_,
495 faceCentres,
496 faceAreas,
497 isChangedCell.toc(),
498 mesh_.cells(),
499 cellCentres,
500 cellVolumes
501 );
502
503 return true;
504 }
505}
506
507
508// ************************************************************************* //
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.
Definition DynamicList.H:68
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
Definition PtrListI.H:289
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
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...
Definition bitSet.H:61
labelList toc() const
The indices of the on bits as a sorted labelList.
Definition bitSet.C:476
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
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.
Definition fvMesh.H:85
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.
Definition polyPatch.H:73
static void updateCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, const List< label > &cellIDs, const List< cell > &cells, vectorField &cellCtrs_s, scalarField &cellVols_s)
Update cell centres and volumes for the cells in the set cellIDs.
static void updateFaceCentresAndAreas(const primitiveMesh &mesh, const UList< label > &faceIDs, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Update face centres and areas for the faces in the set faceIDs.
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).
Definition refPtr.H:54
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.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition fvcMeshPhi.C:29
Namespace for OpenFOAM.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
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
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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)
Definition errorManip.H:139
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).
Definition zero.H:127
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.
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.