Loading...
Searching...
No Matches
displacementPointSmoothingMotionSolver.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) 2024-2025 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 "syncTools.H"
32#include "motionSmootherAlgo.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
39
41 (
45 );
46
48 (
51 displacement
52 );
53}
54
55
56// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57
59(
60 const labelHashSet& changedFaces,
61 labelHashSet& affectedFaces
62)
63{
64 bitSet affectedPoints(mesh().nPoints());
65
66 for (const label facei : changedFaces)
67 {
68 affectedPoints.set(mesh().faces()[facei]);
69 }
70
72 (
73 mesh(),
74 affectedPoints,
75 orEqOp<unsigned int>(),
76 0U
77 );
78
79 for (const label pointi : affectedPoints)
80 {
81 for (const label celli : mesh().pointCells()[pointi])
82 {
83 affectedFaces.insert(mesh().cells()[celli]);
84 }
85 }
86}
87
88
90{
91 if
92 (
93 relaxationFactors_.empty()
94 || (relaxationFactors_.size() == 1 && relaxationFactors_[0] == 1.0)
95 )
96 {
97 relaxedPoints_ = points0() + pointDisplacement().internalField();
98 return true;
99 }
100
101
102 const pointField oldRelaxedPoints(relaxedPoints_);
103
104 labelHashSet affectedFaces(facesToMove_);
105
106 // Create a list of relaxation levels
107 // -1 indicates a point which is not to be moved
108 // 0 is the starting value for a moving point
109 labelList relaxationLevel(mesh().nPoints(), -1);
110
111 for (const label facei : affectedFaces)
112 {
113 for (const label pointi : mesh().faces()[facei])
114 {
115 relaxationLevel[pointi] = 0;
116 }
117 }
118
120 (
121 mesh(),
122 relaxationLevel,
124 label(-1)
125 );
126
127 // Loop whilst relaxation levels are being incremented
128 bool complete(false);
129 while (!complete)
130 {
131 //scalar nAffectedFaces(affectedFaces.size());
132 //reduce(nAffectedFaces, sumOp<scalar>());
133 //Info << " Moving " << nAffectedFaces << " faces" << endl;
134
135 // Move the points
136 forAll(relaxationLevel, pointI)
137 {
138 if (relaxationLevel[pointI] >= 0)
139 {
140 const scalar x
141 (
142 relaxationFactors_[relaxationLevel[pointI]]
143 );
144
145 relaxedPoints_[pointI] =
146 (1 - x)*oldRelaxedPoints[pointI]
147 + x*(points0()[pointI] + pointDisplacement()[pointI]);
148 }
149 }
150
151 // Get a list of changed faces
152 labelHashSet markedFaces;
153 markAffectedFaces(affectedFaces, markedFaces);
154 labelList markedFacesList(markedFaces.toc());
155
156 // Update the geometry
157 meshGeometry_.correct(relaxedPoints_, markedFacesList);
158
159 // Check the modified face quality
160 markedFaces.clear();
162 (
163 false,
164 meshQualityDict_,
165 meshGeometry_,
166 relaxedPoints_,
167 markedFacesList,
168 markedFaces
169 );
170
171 // Mark the affected faces
172 affectedFaces.clear();
173 markAffectedFaces(markedFaces, affectedFaces);
174
175 // Increase relaxation and check convergence
176 bitSet pointsToRelax(mesh().nPoints());
177 complete = true;
178 for (const label facei : affectedFaces)
179 {
180 pointsToRelax.set(mesh().faces()[facei]);
181 }
182
183 for (const label pointi : pointsToRelax)
184 {
185 if (relaxationLevel[pointi] < relaxationFactors_.size() - 1)
186 {
187 ++ relaxationLevel[pointi];
188 complete = false;
189 }
190 }
191
192 // Synchronise relaxation levels
194 (
195 mesh(),
196 relaxationLevel,
198 label(0)
199 );
200
201 // Synchronise completion
202 UPstream::reduceAnd(complete);
203 }
204
205 // Check for convergence
206 bool converged(true);
207 forAll(mesh().faces(), facei)
208 {
209 const face& fPoints(mesh().faces()[facei]);
210
211 for (const label pointi : fPoints)
212 {
213 if (relaxationLevel[pointi] > 0)
214 {
215 facesToMove_.insert(facei);
216
217 converged = false;
218
219 break;
220 }
221 }
222 }
223
224 // Syncronise convergence
225 UPstream::reduceAnd(converged);
226
227 //if (converged)
228 //{
229 // Info<< "... Converged" << endl << endl;
230 //}
231 //else
232 //{
233 // Info<< "... Not converged" << endl << endl;
234 //}
235
236 return converged;
237}
238
239
241(
242 const dictionary& dict
243)
244{
245 if (dict.getOrDefault<bool>("moveInternalFaces", true))
246 {
247 facesToMove_.resize(2*mesh().nFaces());
248 forAll(mesh().faces(), faceI)
249 {
250 facesToMove_.insert(faceI);
251 }
252 }
253 else
254 {
255 facesToMove_.resize(2*(mesh().nBoundaryFaces()));
256 for
257 (
258 label faceI = mesh().nInternalFaces();
259 faceI < mesh().nFaces();
260 ++ faceI
261 )
262 {
263 facesToMove_.insert(faceI);
265 }
266}
267
268
269// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
270
273(
274 const polyMesh& mesh,
275 const IOdictionary& dict
276)
277:
279 meshGeometry_(mesh),
280 pointSmoother_(pointSmoother::New(mesh, coeffDict())),
281 nPointSmootherIter_
282 (
283 coeffDict().get<label>("nPointSmootherIter")
284 ),
285 relaxedPoints_(mesh.points())
286{
287 if (coeffDict().readIfPresent("relaxationFactors", relaxationFactors_))
289 meshQualityDict_ = coeffDict().subDict("meshQuality");
290 }
292}
293
294
297(
298 const polyMesh& mesh,
299 const IOdictionary& dict,
300 const pointVectorField& pointDisplacement,
301 const pointIOField& points0
302)
303:
304 displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
305 meshGeometry_(mesh),
306 pointSmoother_
307 (
308 pointSmoother::New
309 (
310 mesh,
311 coeffDict()
312 )
313 ),
314 nPointSmootherIter_
315 (
316 coeffDict().get<label>("nPointSmootherIter")
317 ),
318 relaxedPoints_(mesh.points())
319{
320 if (coeffDict().readIfPresent("relaxationFactors", relaxationFactors_))
321 {
322 meshQualityDict_ = coeffDict().subDict("meshQuality");
323 }
325}
326
327
328// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
329
333 //Note: twoDCorrect already done by ::solve
334
335 return relaxedPoints_;
336}
337
338
340{
341 //Pout<< "time:" << mesh().time().timeName()
342 // << " mesh faceCentres:" << gAverage(mesh().faceCentres())
343 // << " mesh cellCentres:" << gAverage(mesh().cellCentres())
344 // << endl;
345
346 movePoints(curPoints());
347
348 // Update values on pointDisplacement
349 pointDisplacement().boundaryFieldRef().updateCoeffs();
350
351 // Extend: face-to-point-to-cell-to-faces
352 labelHashSet affectedFaces;
353 markAffectedFaces(facesToMove_, affectedFaces);
354
355 for(label i = 0; i < nPointSmootherIter_; i ++)
356 {
357 meshGeometry_.correct
358 (
359 points0() + pointDisplacement().internalField(),
360 affectedFaces.toc()
361 );
362 //Pout<< "iter:" << i
363 // << " faceCentres:" << gAverage(meshGeometry_.faceCentres())
364 // << " cellCentres:" << gAverage(meshGeometry_.cellCentres())
365 // << endl;
366
367 pointSmoother_->update
368 (
369 affectedFaces.toc(),
370 points0(),
371 points0() + pointDisplacement().internalField(),
372 meshGeometry_,
373 pointDisplacement()
374 );
375 }
376
377 relax();
378
379 twoDCorrectPoints(relaxedPoints_);
380
381 // Update pointDisplacement for actual relaxedPoints. Keep fixed-value
382 // bcs.
383 pointDisplacement().primitiveFieldRef() = relaxedPoints_-points0();
384
385 // Adhere to multi-point constraints
386 const pointConstraints& pcs =
387 pointConstraints::New(pointDisplacement().mesh());
388 pcs.constrainDisplacement(pointDisplacement(), false);
389
390 // Update relaxedPoints to take constraints into account
391 relaxedPoints_ = points0() + pointDisplacement().internalField();
392}
393
394
395// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
propsDict readIfPresent("fields", acceptFields)
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
void clear()
Remove all entries from table.
Definition HashTable.C:742
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
static void reduceAnd(bool &value, const int communicator=worldComm)
Logical (and) reduction (MPI_AllReduce).
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
Virtual base class for displacement motion solver.
pointVectorField & pointDisplacement() noexcept
Return reference to the point motion displacement field.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
Quality-based under-relaxation for run-time selectable point smoothing.
autoPtr< pointSmoother > pointSmoother_
Point smoothing method.
displacementPointSmoothingMotionSolver(const polyMesh &, const IOdictionary &)
Construct from a polyMesh and an IOdictionary.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
labelHashSet facesToMove_
Set of the faces which are to be moved.
polyMeshGeometry meshGeometry_
Part-updatable mesh geometry.
scalarList relaxationFactors_
Relaxation factors to use in each iteration.
void markAffectedFaces(const labelHashSet &changedFaces, labelHashSet &affectedFaces)
Mark affected faces.
const label nPointSmootherIter_
Number of point smoother iterations per timestep.
virtual void setFacesToMove(const dictionary &)
Set all the faces to be moved.
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
virtual void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition pointCells.H:55
Application of (multi-)patch point constraints.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Abstract base class for point smoothing methods. Handles parallel communication via reset and average...
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
label nFaces() const noexcept
Number of mesh faces.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
UEqn relax()
dynamicFvMesh & mesh
const pointField & points
label nPoints
const cellShapeList & cells
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Namespace for OpenFOAM.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
vectorField pointField
pointField is a vectorField.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))