Loading...
Searching...
No Matches
displacementSBRStressFvMotionSolver.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2017 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
30#include "motionInterpolation.H"
31#include "motionDiffusivity.H"
32#include "fvmLaplacian.H"
34#include "fvcDiv.H"
35#include "fvcGrad.H"
36#include "surfaceInterpolate.H"
37#include "fvcLaplacian.H"
38#include "mapPolyMesh.H"
39#include "fvOptions.H"
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
46
48 (
52 );
53
55 (
58 displacement
59 );
60}
61
62
63// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64
65Foam::displacementSBRStressFvMotionSolver::displacementSBRStressFvMotionSolver
66(
67 const polyMesh& mesh,
68 const IOdictionary& dict
69)
70:
73 cellDisplacement_
74 (
76 (
77 "cellDisplacement",
78 mesh.time().timeName(),
79 mesh,
80 IOobject::READ_IF_PRESENT,
81 IOobject::AUTO_WRITE
82 ),
83 fvMesh_,
84 dimensionedVector(pointDisplacement().dimensions(), Zero),
85 cellMotionBoundaryTypes<vector>(pointDisplacement().boundaryField())
86 ),
87 interpolationPtr_
88 (
89 coeffDict().found("interpolation")
90 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
91 : motionInterpolation::New(fvMesh_)
92 ),
93 diffusivityPtr_
94 (
95 motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
96 )
97{}
98
99
100Foam::displacementSBRStressFvMotionSolver::
101displacementSBRStressFvMotionSolver
102(
103 const polyMesh& mesh,
104 const IOdictionary& dict,
105 const pointVectorField& pointDisplacement,
106 const pointIOField& points0
107)
108:
109 displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
111 cellDisplacement_
112 (
114 (
115 "cellDisplacement",
116 mesh.time().timeName(),
117 mesh,
118 IOobject::READ_IF_PRESENT,
119 IOobject::AUTO_WRITE
120 ),
121 fvMesh_,
123 (
124 displacementMotionSolver::pointDisplacement().dimensions(),
125 Zero
126 ),
127 cellMotionBoundaryTypes<vector>
128 (
129 displacementMotionSolver::pointDisplacement().boundaryField()
130 )
131 ),
132 interpolationPtr_
133 (
134 coeffDict().found("interpolation")
135 ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
136 : motionInterpolation::New(fvMesh_)
137 ),
138 diffusivityPtr_
139 (
141 )
142{}
143
144
145// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146
149{}
150
151
152// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153
156{
157 interpolationPtr_->interpolate
158 (
159 cellDisplacement_,
160 pointDisplacement_
161 );
162
163 // Evaluate the bcs so they are consistent with the internal field
164 // Might fight the multi-patch behaviour inside volPointInterpolate
165 if
166 (
167 pointDisplacement_.boundaryField().size()
168 != cellDisplacement_.boundaryField().size()
169 )
170 {
171 pointDisplacement_.correctBoundaryConditions();
172 }
173
174 tmp<pointField> tcurPoints
175 (
176 points0() + pointDisplacement().primitiveField()
177 );
179 twoDCorrectPoints(tcurPoints.ref());
180
181 return tcurPoints;
182}
183
184
186{
187 // The points have moved so before interpolation update
188 // the motionSolver accordingly
189 movePoints(fvMesh_.points());
190
191 diffusivityPtr_->correct();
192 pointDisplacement_.boundaryFieldRef().updateCoeffs();
193
194 const surfaceScalarField Df
195 (
196 dimensionedScalar("viscosity", dimViscosity, 1.0)
197 *diffusivityPtr_->operator()()
198 );
199
200 // Make sure the cellMotion bcs are consistent with the pointDisplacement.
201 // This makes sure the grad below uses more up-to-date values.
202 cellDisplacement_.boundaryFieldRef().updateCoeffs();
203
204 volTensorField gradCd("gradCd", fvc::grad(cellDisplacement_));
205
207
209 (
211 (
212 2*Df,
213 cellDisplacement_,
214 "laplacian(diffusivity,cellDisplacement)"
215 )
216
217 + fvc::div
218 (
219 Df
220 *(
222 (
223 cellDisplacement_.mesh().Sf(),
224 gradCd.T() - gradCd
225 )
226
227 // Solid-body rotation "lambda" term
228 - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
229 )
230 )
231
232 /*
233 - fvc::laplacian
234 (
235 2*Df,
236 cellDisplacement_,
237 "laplacian(diffusivity,cellDisplacement)"
238 )
239
240 + fvc::div
241 (
242 Df
243 *(
244 fvc::dotInterpolate
245 (
246 cellDisplacement_.mesh().Sf(),
247 gradCd + gradCd.T()
248 )
249
250 // Solid-body rotation "lambda" term
251 - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
252 )
253 )
254 */
255 ==
256 fvOptions(cellDisplacement_)
257 );
259 fvOptions.constrain(TEqn);
260 TEqn.solveSegregatedOrCoupled();
261 fvOptions.correct(cellDisplacement_);
262}
263
264
266(
267 const mapPolyMesh& mpm
268)
269{
271
272 // Update diffusivity. Note two stage to make sure old one is de-registered
273 // before creating/registering new one.
274 diffusivityPtr_.reset(nullptr);
275 diffusivityPtr_ = motionDiffusivity::New
276 (
277 fvMesh_,
278 coeffDict().lookup("diffusivity")
279 );
280}
281
282
283// ************************************************************************* //
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
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field).
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ 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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
Mesh motion solver for an fvMesh. Based on solving the cell-centre solid-body rotation stress equatio...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void updateMesh(const mapPolyMesh &)
Update topology.
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 void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
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
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the laplacian of the given field.
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
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
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.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dictionary dict
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))