displacementMotionSolverMeshMover.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2015 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "pointConstraints.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38
40 (
44 );
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
51(
52 const dictionary& moveDict,
53 const label nAllowableErrors,
54 labelList& checkFaces
55)
56{
57 const label nRelaxIter = moveDict.get<label>("nRelaxIter");
58
59 meshMover_.setDisplacementPatchFields();
60
61 Info<< typeName << " : Moving mesh ..." << endl;
62
63 scalar oldErrorReduction = -1;
64
65 bool meshOk = false;
66
67 for (label iter = 0; iter < 2*nRelaxIter; ++ iter)
68 {
69 Info<< typeName << " : Iteration " << iter << endl;
70
71 if (iter == nRelaxIter)
72 {
73 Info<< typeName
74 << " : Displacement scaling for error reduction set to 0."
75 << endl;
76 oldErrorReduction = meshMover_.setErrorReduction(0.0);
77 }
78
79 if
80 (
81 meshMover_.scaleMesh
82 (
83 checkFaces,
85 meshMover_.paramDict(),
86 moveDict,
87 true,
88 nAllowableErrors
89 )
90 )
91 {
92 Info<< typeName << " : Successfully moved mesh" << endl;
93 meshOk = true;
94 break;
95 }
96 }
97
98 if (oldErrorReduction >= 0)
99 {
100 meshMover_.setErrorReduction(oldErrorReduction);
101 }
102
103 Info<< typeName << " : Finished moving mesh ..." << endl;
104
105 return meshOk;
106}
107
108
109// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110
112(
113 const dictionary& dict,
114 const List<labelPair>& baffles,
115 pointVectorField& pointDisplacement,
116 const bool dryRun
117)
118:
119 externalDisplacementMeshMover(dict, baffles, pointDisplacement, dryRun),
120
121 solverPtr_
122 (
124 (
125 dict.get<word>("solver"),
126 pointDisplacement.mesh()(),
128 (
130 (
131 "motionSolverDict",
132 pointDisplacement.mesh().time().constant(),
133 pointDisplacement.db(),
134 IOobject::NO_READ,
135 IOobject::NO_WRITE,
136 false
137 ),
138 dict
139 ),
140 pointDisplacement,
142 (
144 (
145 "points0",
146 pointDisplacement.mesh().time().constant(),
147 pointDisplacement.db(),
148 IOobject::NO_READ,
149 IOobject::NO_WRITE,
150 false
151 ),
152 pointDisplacement.mesh()().points()
153 )
154 )
155 ),
156
157 adaptPatchIDs_(getFixedValueBCs(pointDisplacement)),
158 adaptPatchPtr_(getPatch(mesh(), adaptPatchIDs_)),
159
160 scale_
161 (
163 (
164 "scale",
165 pointDisplacement.time().timeName(),
166 pointDisplacement.db(),
167 IOobject::NO_READ,
168 IOobject::AUTO_WRITE
169 ),
170 pMesh(),
171 dimensionedScalar("scale", dimless, 1.0)
172 ),
173
174 oldPoints_(mesh().points()),
175
176 meshMover_
177 (
178 const_cast<polyMesh&>(mesh()),
179 const_cast<pointMesh&>(pMesh()),
180 adaptPatchPtr_(),
181 pointDisplacement,
182 scale_,
183 oldPoints_,
184 adaptPatchIDs_,
185 dict,
186 dryRun
187 ),
188
189 fieldSmoother_(mesh())
190{}
191
192
193// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
194
196{}
197
198
199// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
200
202(
203 const dictionary& moveDict,
204 const label nAllowableErrors,
205 labelList& checkFaces
206)
207{
208 // Correct and smooth the patch displacements so points next to
209 // points where the extrusion was disabled also use less extrusion.
210 // Note that this has to update the pointDisplacement boundary conditions
211 // as well, not just the internal field.
212 {
213 const label nSmoothPatchThickness = meshRefinement::get<label>
214 (
215 moveDict, "nSmoothThickness", dryRun_, keyType::REGEX
216 );
217
218 const word minThicknessName = meshRefinement::get<word>
219 (
220 moveDict, "minThicknessName", dryRun_, keyType::REGEX, word::null
221 );
222
223 scalarField zeroMinThickness;
224
225 if (minThicknessName == "none")
226 {
227 zeroMinThickness = scalarField(adaptPatchPtr_().nPoints(), Zero);
228 }
229
230 const scalarField& minThickness =
231 (
232 (minThicknessName == "none")
233 ? zeroMinThickness
234 : mesh().lookupObject<scalarField>(minThicknessName)
235 );
236
237 const bitSet isPatchMasterPoint
238 (
240 (
241 mesh(),
242 adaptPatchPtr_().meshPoints()
243 )
244 );
245
246 const bitSet isPatchMasterEdge
247 (
249 (
250 mesh(),
251 adaptPatchPtr_().meshEdges
252 (
253 mesh().edges(),
254 mesh().pointEdges()
255 )
256 )
257 );
258
259 // Smooth patch displacement
260
261 vectorField displacement
262 (
263 pointDisplacement().internalField(),
264 adaptPatchPtr_().meshPoints()
265 );
266
267 fieldSmoother_.minSmoothField
268 (
269 nSmoothPatchThickness,
270 isPatchMasterPoint,
271 isPatchMasterEdge,
272 adaptPatchPtr_(),
273 minThickness,
274 displacement
275 );
276
277
278 scalar resid = 0;
279
280 forAll(displacement, patchPointI)
281 {
282 const label pointI(adaptPatchPtr_().meshPoints()[patchPointI]);
283
284 resid += mag(pointDisplacement()[pointI]-displacement[patchPointI]);
285
286 pointDisplacement()[pointI] = displacement[patchPointI];
287 }
288
289 // Take over smoothed displacements on bcs
290 meshMover_.setDisplacementPatchFields();
291 }
292
293 // Use motionSolver to calculate internal displacement
294 {
295 solverPtr_->pointDisplacement() == pointDisplacement();
296 // Force solving and constraining - just so its pointDisplacement gets
297 // the correct value
298 (void)solverPtr_->newPoints();
299 pointDisplacement() == solverPtr_->pointDisplacement();
300 }
301
302 return moveMesh(moveDict, nAllowableErrors, checkFaces);
303}
304
305
307{
309
310 // Update motion solver for new geometry
311 solverPtr_->movePoints(p);
312
313 // Update motionSmoother for new geometry (moves adaptPatchPtr_)
314 meshMover_.movePoints();
315
316 // Assume current mesh location is correct (reset oldPoints, scale)
317 meshMover_.correct();
318}
319
320
321// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Quality-based under-relaxation wrapped around generic displacementMotionSolver.
Virtual base class for displacement motion solver.
virtual void move()=0
Virtual base class for mesh movers with externally provided displacement field giving the boundary co...
List< labelPair > baffles_
Baffles in the mesh.
@ REGEX
Regular expression.
Definition: keyType.H:82
static bitSet getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
static bitSet getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
scalar setErrorReduction(const scalar)
Set the errorReduction (by how much to scale the displacement.
void movePoints()
Update for new mesh geometry.
const dictionary & paramDict() const
static void setDisplacementPatchFields(const labelList &patchIDs, pointVectorField &pointDisplacement)
Set patch fields on patchIDs to be consistent with.
bool scaleMesh(labelList &checkFaces, const bool smoothMesh=true, const label nAllow=0)
Move mesh with given scale. Return true if mesh ok or has.
const Type & lookupObject(const word &name, const bool recursive=false) const
virtual void moveMesh()
constant condensation/saturation model.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
const pointField & points
label nPoints
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333