solidBodyDisplacementLaplacianFvMotionSolver.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) 2016-2020 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
29#include "motionInterpolation.H"
30#include "motionDiffusivity.H"
31#include "fvmLaplacian.H"
33#include "OFstream.H"
34#include "meshTools.H"
35#include "mapPolyMesh.H"
37#include "transformField.H"
38#include "fvOptions.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
45
47 (
51 );
52
54 (
57 displacement
58 );
59}
60
61
62// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63
64Foam::solidBodyDisplacementLaplacianFvMotionSolver::
65solidBodyDisplacementLaplacianFvMotionSolver
66(
67 const polyMesh& mesh,
68 const IOdictionary& dict
69)
70:
73 SBMFPtr_(solidBodyMotionFunction::New(coeffDict(), mesh.time())),
74 cellDisplacement_
75 (
77 (
78 "cellDisplacement",
79 mesh.time().timeName(),
80 mesh,
81 IOobject::READ_IF_PRESENT,
82 IOobject::AUTO_WRITE
83 ),
84 fvMesh_,
85 dimensionedVector(pointDisplacement_.dimensions(), Zero),
86 cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
87 ),
88 pointLocation_(nullptr),
89 interpolationPtr_
90 (
91 coeffDict().found("interpolation")
92 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
93 : motionInterpolation::New(fvMesh_)
94 ),
95 diffusivityPtr_
96 (
97 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
98 ),
99 frozenPointsZone_
100 (
101 coeffDict().found("frozenPointsZone")
102 ? fvMesh_.pointZones().findZoneID
103 (
104 coeffDict().get<word>("frozenPointsZone")
105 )
106 : -1
107 )
108{
110 (
111 "pointLocation",
113 fvMesh_,
116 );
117
118 if (debug)
119 {
120 Info<< "solidBodyDisplacementLaplacianFvMotionSolver:" << nl
121 << " diffusivity : " << diffusivityPtr_().type() << nl
122 << " frozenPoints zone : " << frozenPointsZone_ << endl;
123 }
124
125
127 {
128 pointLocation_.reset
129 (
131 (
132 io,
134 )
135 );
136
137 if (debug)
138 {
139 Info<< "solidBodyDisplacementLaplacianFvMotionSolver :"
140 << " Read pointVectorField "
141 << io.name()
142 << " to be used for boundary conditions on points."
143 << nl
144 << "Boundary conditions:"
145 << pointLocation_().boundaryField().types() << endl;
146 }
147 }
148}
149
150
151Foam::solidBodyDisplacementLaplacianFvMotionSolver::
152solidBodyDisplacementLaplacianFvMotionSolver
153(
154 const polyMesh& mesh,
155 const IOdictionary& dict,
156 const pointVectorField& pointDisplacement,
157 const pointIOField& points0
158)
159:
160 displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
162 SBMFPtr_(solidBodyMotionFunction::New(coeffDict(), mesh.time())),
163 cellDisplacement_
164 (
166 (
167 "cellDisplacement",
168 mesh.time().timeName(),
169 mesh,
170 IOobject::READ_IF_PRESENT,
171 IOobject::AUTO_WRITE
172 ),
173 fvMesh_,
174 dimensionedVector(pointDisplacement_.dimensions(), Zero),
175 cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
176 ),
177 pointLocation_(nullptr),
178 interpolationPtr_
179 (
180 coeffDict().found("interpolation")
181 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
182 : motionInterpolation::New(fvMesh_)
183 ),
184 diffusivityPtr_
185 (
186 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
187 ),
188 frozenPointsZone_
189 (
190 coeffDict().found("frozenPointsZone")
191 ? fvMesh_.pointZones().findZoneID
192 (
193 coeffDict().get<word>("frozenPointsZone")
194 )
195 : -1
196 )
197{
199 (
200 "pointLocation",
202 fvMesh_,
205 );
206
207 if (debug)
208 {
209 Info<< "solidBodyDisplacementLaplacianFvMotionSolver:" << nl
210 << " diffusivity : " << diffusivityPtr_().type() << nl
211 << " frozenPoints zone : " << frozenPointsZone_ << endl;
212 }
213
214
216 {
217 pointLocation_.reset
218 (
220 (
221 io,
223 )
224 );
225
226 if (debug)
227 {
228 Info<< "solidBodyDisplacementLaplacianFvMotionSolver :"
229 << " Read pointVectorField "
230 << io.name()
231 << " to be used for boundary conditions on points."
232 << nl
233 << "Boundary conditions:"
234 << pointLocation_().boundaryField().types() << endl;
235 }
236 }
237}
238
239
240// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
241
244{}
245
246
247// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
248
251{
252 if (!diffusivityPtr_)
253 {
254 diffusivityPtr_ = motionDiffusivity::New
255 (
256 fvMesh_,
257 coeffDict().lookup("diffusivity")
258 );
259 }
260
261 return *diffusivityPtr_;
262}
263
264
267{
268 interpolationPtr_->interpolate
269 (
270 cellDisplacement_,
271 pointDisplacement_
272 );
273
274 tmp<pointField> tnewPoints
275 (
276 transformPoints(SBMFPtr_().transformation(), points0())
277 );
278 const pointField& newPoints = tnewPoints();
279
280 if (pointLocation_)
281 {
282 if (debug)
283 {
284 Info<< "solidBodyDisplacementLaplacianFvMotionSolver : applying "
285 << " boundary conditions on " << pointLocation_().name()
286 << " to new point location."
287 << endl;
288 }
289
290 pointLocation_().primitiveFieldRef() =
291 newPoints
292 + pointDisplacement_.internalField();
293
294 pointLocation_().correctBoundaryConditions();
295
296 // Implement frozen points
297 if (frozenPointsZone_ != -1)
298 {
299 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
300
301 forAll(pz, i)
302 {
303 pointLocation_()[pz[i]] = newPoints[pz[i]];
304 }
305 }
306
307 twoDCorrectPoints(pointLocation_().primitiveFieldRef());
308
309 return tmp<pointField>(pointLocation_().primitiveField());
310 }
311 else
312 {
313 tmp<pointField> tcurPoints
314 (
315 newPoints + pointDisplacement_.primitiveField()
316 );
317 pointField& curPoints = tcurPoints.ref();
318
319 // Implement frozen points
320 if (frozenPointsZone_ != -1)
321 {
322 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
323
324 forAll(pz, i)
325 {
326 curPoints[pz[i]] = newPoints[pz[i]];
327 }
328 }
329
330 twoDCorrectPoints(curPoints);
331
332 return tcurPoints;
333 }
334}
335
336
338{
339 // The points have moved so before interpolation update
340 // the motionSolver accordingly
341 movePoints(fvMesh_.points());
342
343 diffusivity().correct();
344 pointDisplacement_.boundaryFieldRef().updateCoeffs();
345
347
349 (
351 (
352 dimensionedScalar("viscosity", dimViscosity, 1.0)
353 *diffusivity().operator()(),
354 cellDisplacement_,
355 "laplacian(diffusivity,cellDisplacement)"
356 )
357 ==
358 fvOptions(cellDisplacement_)
359 );
360
363 fvOptions.correct(cellDisplacement_);
364}
365
366
368(
369 const mapPolyMesh& mpm
370)
371{
373
374 // Update diffusivity. Note two stage to make sure old one is de-registered
375 // before creating/registering new one.
376 diffusivityPtr_.clear();
377}
378
379
380// ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
fv::options & fvOptions
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
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Virtual base class for displacement motion solver.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:62
const dictionary & solverDict() const
Return the solver dictionary taking into account finalIteration.
Definition: fvMatrix.C:1558
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Finite-volume options.
Definition: fvOptions.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Abstract base class for cell-centre mesh motion diffusivity.
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
void updateMesh()
Update for new mesh topology.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:61
A subset of mesh points.
Definition: pointZone.H:68
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Lookup type of boundary radiation properties.
Definition: lookup.H:66
Applies Laplacian displacement solving on top of a transformation of the initial points using a solid...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
motionDiffusivity & diffusivity()
Return reference to the diffusivity field.
Base class for defining solid-body motions.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
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
dynamicFvMesh & mesh
Calculate the matrix for the laplacian of the field.
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Namespace for OpenFOAM.
const dimensionSet dimViscosity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
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.
void transformPoints(vectorField &, const septernion &, const vectorField &)
Transform given vectorField of coordinates with the given septernion.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Spatial transformation functions for primitive fields.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
conserve primitiveFieldRef()+