Loading...
Searching...
No Matches
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"
38#include "fvOptions.H"
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
126 if (io.typeHeaderOk<pointVectorField>(true))
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{
198 IOobject io
199 (
200 "pointLocation",
201 fvMesh_.time().timeName(),
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
215 if (io.typeHeaderOk<pointVectorField>(true))
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;
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 );
260
261 return *diffusivityPtr_;
262}
263
264
267{
268 interpolationPtr_->interpolate
269 (
270 cellDisplacement_,
271 pointDisplacement_
272 );
273
274 // Evaluate the bcs so they are consistent with the internal field
275 // Might fight the multi-patch behaviour inside volPointInterpolate
276 if
277 (
278 pointDisplacement_.boundaryField().size()
279 != cellDisplacement_.boundaryField().size()
280 )
281 {
282 pointDisplacement_.correctBoundaryConditions();
283 }
284
285 tmp<pointField> tnewPoints
286 (
287 transformPoints(SBMFPtr_().transformation(), points0())
288 );
289 const pointField& newPoints = tnewPoints();
290
291 if (pointLocation_)
292 {
293 if (debug)
294 {
295 Info<< "solidBodyDisplacementLaplacianFvMotionSolver : applying "
296 << " boundary conditions on " << pointLocation_().name()
297 << " to new point location."
298 << endl;
299 }
300
301 pointLocation_().primitiveFieldRef() =
302 newPoints
303 + pointDisplacement_.internalField();
304
305 pointLocation_().correctBoundaryConditions();
306
307 // Implement frozen points
308 if (frozenPointsZone_ != -1)
309 {
310 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
311
312 forAll(pz, i)
313 {
314 pointLocation_()[pz[i]] = newPoints[pz[i]];
315 }
316 }
317
318 twoDCorrectPoints(pointLocation_().primitiveFieldRef());
319
320 return tmp<pointField>(pointLocation_().primitiveField());
321 }
322 else
323 {
324 tmp<pointField> tcurPoints
325 (
326 newPoints + pointDisplacement_.primitiveField()
327 );
328 pointField& curPoints = tcurPoints.ref();
329
330 // Implement frozen points
331 if (frozenPointsZone_ != -1)
332 {
333 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
334
335 forAll(pz, i)
336 {
337 curPoints[pz[i]] = newPoints[pz[i]];
338 }
339 }
340
342
343 return tcurPoints;
344 }
345}
346
347
349{
350 // The points have moved so before interpolation update
351 // the motionSolver accordingly
352 movePoints(fvMesh_.points());
353
354 diffusivity().correct();
355 pointDisplacement_.boundaryFieldRef().updateCoeffs();
356
357 fv::options& fvOptions(fv::options::New(fvMesh_));
358
360 (
362 (
363 dimensionedScalar("viscosity", dimViscosity, 1.0)
364 *diffusivity().operator()(),
365 cellDisplacement_,
366 "laplacian(diffusivity,cellDisplacement)"
367 )
368 ==
369 fvOptions(cellDisplacement_)
370 );
372 fvOptions.constrain(TEqn);
373 TEqn.solveSegregatedOrCoupled();
374 fvOptions.correct(cellDisplacement_);
375}
376
377
379(
380 const mapPolyMesh& mpm
381)
382{
384
385 // Update diffusivity. Note two stage to make sure old one is de-registered
386 // before creating/registering new one.
387 diffusivityPtr_.clear();
388}
389
390
391// ************************************************************************* //
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
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
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.
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.
pointVectorField pointDisplacement_
Point motion field.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
wordList cellMotionBoundaryTypes(const typename GeometricField< Type, pointPatchField, pointMesh >::Boundary &pmUbf) const
Create the corresponding patch types for cellMotion from those.
fvMotionSolver(const polyMesh &)
Construct from polyMesh.
const fvMesh & mesh() const
Return reference to the fvMesh to be moved.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Abstract base class for cell-centre mesh motion diffusivity.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
Virtual base class for mesh motion solver.
virtual tmp< pointField > newPoints()
Provide new points for motion. Solves for motion.
virtual void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
A subset of mesh points.
Definition pointZone.H:64
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Lookup type of boundary radiation properties.
Definition lookup.H:60
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:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
const auto & io
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))
word timeName
Definition getTimeIndex.H:3
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 handling debugging switches.
Definition debug.C:45
int debug
Static debugging option.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Namespace for OpenFOAM.
const dimensionSet dimViscosity
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.
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
void transformPoints(vectorField &, const septernion &, const vectorField &)
Transform given vectorField of coordinates with the given septernion.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Spatial transformation functions for primitive fields.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
conserve primitiveFieldRef()+