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"
31#include "primitiveMeshTools.H"
32#include "emptyPolyPatch.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
43 dict
44 );
45}
46
47
48void Foam::solidBodyFvGeometryScheme::setMeshMotionData()
49{
50 if (!cacheInitialised_ || !cacheMotion_)
51 {
52 DebugInFunction << "Creating cache" << endl;
53
54 changedFaceIDs_.clear(); // used for face areas, meshPhi
55 changedPatchIDs_.clear(); // used for meshPhi
56 changedCellIDs_.clear(); // used for cell volumes
57
58 const pointField& oldPoints = mesh_.oldPoints();
59 const pointField& currPoints = mesh_.points();
60
61 if (oldPoints.size() != currPoints.size())
62 {
64 << "Old and current points sizes must be the same. "
65 << "Old points:" << oldPoints.size()
66 << " Current points:" << currPoints.size()
67 << abort(FatalError);
68 }
69
70 bitSet changedPoints(oldPoints.size());
71
72 // Check for non-identical points
73 forAll(changedPoints, pointi)
74 {
75 changedPoints.set(pointi, oldPoints[pointi] != currPoints[pointi]);
76 }
77
79 << "SBM --- Changed points:"
80 << returnReduce(changedPoints.count(), sumOp<label>())
81 << endl;
82
83 // Quick return if no points have moved
84 if (returnReduce(changedPoints.none(), andOp<label>()))
85 {
86 return;
87 }
88
89 bitSet cellIDs(mesh_.nCells());
90 bitSet faceIDs(mesh_.nFaces());
91
92 const auto& pointFaces = mesh_.pointFaces();
93 const auto& own = mesh_.faceOwner();
94 const auto& nbr = mesh_.faceNeighbour();
95
96 // Identify faces and cells attached to moving points
97 for (const label pointi : changedPoints)
98 {
99 for (const auto facei : pointFaces[pointi])
100 {
101 faceIDs.set(facei);
102
103 cellIDs.set(own[facei]);
104 if (facei < mesh_.nInternalFaces())
105 {
106 cellIDs.set(nbr[facei]);
107 }
108 }
109 }
110
111 changedCellIDs_ = cellIDs.toc();
112
114 << "SBM --- Changed cells:"
115 << returnReduce(changedCellIDs_.size(), sumOp<label>())
116 << endl;
117
118
119 // Construct face and patch ID info
120
121 const auto changedFaceFlag = faceIDs.values();
122
123 DynamicList<label> changedFaceIDs(faceIDs.count());
124 DynamicList<label> changedPatchIDs(faceIDs.count());
125 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
126 {
127 if (changedFaceFlag[facei])
128 {
129 changedFaceIDs.append(facei);
130 changedPatchIDs.append(-1);
131 }
132 }
133
134 const auto& pbm = mesh_.boundaryMesh();
135 for (label patchi = 0; patchi < pbm.size(); ++patchi)
136 {
137 const polyPatch& pp = pbm[patchi];
138
139 for (const label meshFacei : pp.range())
140 {
141 if (changedFaceFlag[meshFacei])
142 {
143 changedFaceIDs.append(meshFacei);
144 changedPatchIDs.append(patchi);
145 }
146 }
147 }
148
149 changedFaceIDs_.transfer(changedFaceIDs);
150 changedPatchIDs_.transfer(changedPatchIDs);
151 }
152
153 cacheInitialised_ = true;
154}
155
156
157// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
158
160(
161 const fvMesh& mesh,
162 const dictionary& dict
163)
164:
166 partialUpdate_(dict.getOrDefault<bool>("partialUpdate", true)),
167 cacheMotion_(dict.getOrDefault<bool>("cacheMotion", true)),
168 cacheInitialised_(false),
169 changedFaceIDs_(),
170 changedPatchIDs_(),
171 changedCellIDs_()
172{
174 << "partialUpdate:" << partialUpdate_
175 << " cacheMotion:" << cacheMotion_
176 << endl;
177}
178
179
180// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181
183{
184 // Note: not calling fvGeometryScheme::movePoints since we want to perform
185 // our own geometry manipulations
186
187 bool haveGeometry =
188 mesh_.hasCellCentres()
189 && mesh_.hasFaceCentres()
190 && mesh_.hasCellVolumes()
191 && mesh_.hasFaceAreas();
192
193 if (!haveGeometry)
194 {
196 << "Creating initial geometry using primitiveMesh::updateGeom"
197 << endl;
198
199 const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
200 return;
201 }
202
203 if (mesh_.moving())
204 {
205 setMeshMotionData();
206
207 DebugInFunction << "Performing partial meshPhi construction" << endl;
208
209 const pointField& oldPoints = mesh_.oldPoints();
210 const pointField& currPoints = mesh_.points();
211
212 if (oldPoints.size() != currPoints.size())
213 {
215 << "Old and current points sizes must be the same. "
216 << "Old points:" << oldPoints.size()
217 << " Current points:" << currPoints.size()
218 << abort(FatalError);
219 }
220
221 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
222 const faceList& faces = mesh_.faces();
223
224 auto tmeshPhi(const_cast<fvMesh&>(mesh_).setPhi());
225 if (tmeshPhi)
226 {
227 const scalar rdt = 1.0/mesh_.time().deltaTValue();
228
229 // Set the mesh flux
230 auto& meshPhi = tmeshPhi.ref();
231 auto& meshPhii = meshPhi.primitiveFieldRef();
232 auto& meshPhiBf = meshPhi.boundaryFieldRef();
233
234 //meshPhi == dimensionedScalar(dimVolume/dimTime, Zero);
235 meshPhii = Zero;
236 meshPhiBf == Zero;
237
238 forAll(changedFaceIDs_, i)
239 {
240 const face& f = faces[changedFaceIDs_[i]];
241
242 if (changedPatchIDs_[i] == -1)
243 {
244 const label facei = changedFaceIDs_[i];
245 meshPhii[facei] = f.sweptVol(oldPoints, currPoints)*rdt;
246 }
247 else
248 {
249 const label patchi = changedPatchIDs_[i];
250 const polyPatch& pp = pbm[patchi];
251
252 if (isA<emptyPolyPatch>(pp))
253 {
254 continue;
255 }
256
257 const label patchFacei = changedFaceIDs_[i] - pp.start();
258
259 meshPhiBf[patchi][patchFacei] =
260 f.sweptVol(oldPoints, currPoints)*rdt;
261 }
262 }
263 }
264
265 if (partialUpdate_ && haveGeometry)
266 {
267 // Keep base geometry and update as needed
268 DebugInFunction << "Performing partial geometry update" << endl;
269
270 // Initialise geometry using the old/existing values
271 vectorField faceCentres(mesh_.faceCentres());
272 vectorField faceAreas(mesh_.faceAreas());
273
274 // Make face centres and areas consistent with new points
276 (
277 mesh_,
278 changedFaceIDs_,
279 mesh_.points(),
280 faceCentres,
281 faceAreas
282 );
283
284 vectorField cellCentres(mesh_.cellCentres());
285 scalarField cellVolumes(mesh_.cellVolumes());
286
288 (
289 mesh_,
290 faceCentres,
291 faceAreas,
292 changedCellIDs_,
293 mesh_.cells(),
294 cellCentres,
295 cellVolumes
296 );
297
298 const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
299 (
300 std::move(faceCentres),
301 std::move(faceAreas),
302 std::move(cellCentres),
303 std::move(cellVolumes)
304 );
305
306 if (debug)
307 {
308 for (const auto& p : mesh_.boundaryMesh())
309 {
310 Pout<< "SBM --- " << p.name()
311 << " sum(Sf)=" << sum(p.faceAreas())
312 << " sum(mag(Sf))=" << sum(mag(p.faceAreas()))
313 << endl;
314 }
315 }
316 }
317 else
318 {
320 << "Performing complete geometry clear and update" << endl;
321
322 // Clear out old geometry
323 // Note: this recreates the old primitiveMesh::movePoints behaviour
324 const_cast<fvMesh&>(mesh_).primitiveMesh::clearGeom();
325
326 // Use lower level to calculate the geometry
327 const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
328 }
329 }
330 else
331 {
332 DebugInFunction << "Performing complete geometry update" << endl;
333
334 // Use lower level to calculate the geometry
335 const_cast<fvMesh&>(mesh_).primitiveMesh::updateGeom();
336 }
337}
338
339
341{
342 cacheInitialised_ = false;
343}
344
345
346// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
labelList cellIDs
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:330
void transfer(List< T > &list)
Definition: List.C:447
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
void append(T *ptr)
Append an element to the end of the list.
Definition: PtrListI.H:113
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Default geometry calculation scheme. Slight stabilisation for bad meshes.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Abstract base class for geometry calculation schemes.
const fvMesh & mesh_
Hold reference to mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
void updateMesh()
Update for new mesh topology.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
virtual const pointField & oldPoints() const
Return old points (mesh motion)
Definition: polyMesh.C:1133
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
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.
label nInternalFaces() const noexcept
Number of internal faces.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const labelListList & pointFaces() const
virtual void updateGeom()
Update all geometric data.
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
Geometry calculation scheme that performs geometry updates only in regions where the mesh has changed...
virtual void movePoints()
Do what is necessary if the mesh has moved.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
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.
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.