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 -------------------------------------------------------------------------------
11 License
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"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(displacementSBRStressFvMotionSolver, 0);
46 
48  (
49  motionSolver,
50  displacementSBRStressFvMotionSolver,
51  dictionary
52  );
53 
55  (
56  displacementMotionSolver,
57  displacementSBRStressFvMotionSolver,
58  displacement
59  );
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
65 Foam::displacementSBRStressFvMotionSolver::displacementSBRStressFvMotionSolver
66 (
67  const polyMesh& mesh,
68  const IOdictionary& dict
69 )
70 :
73  cellDisplacement_
74  (
75  IOobject
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 
100 Foam::displacementSBRStressFvMotionSolver::
101 displacementSBRStressFvMotionSolver
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  (
113  IOobject
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  (
140  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
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  tmp<pointField> tcurPoints
164  (
165  points0() + pointDisplacement().primitiveField()
166  );
168  twoDCorrectPoints(tcurPoints.ref());
169 
170  return tcurPoints;
171 }
172 
173 
175 {
176  // The points have moved so before interpolation update
177  // the motionSolver accordingly
178  movePoints(fvMesh_.points());
179 
180  diffusivityPtr_->correct();
181  pointDisplacement_.boundaryFieldRef().updateCoeffs();
182 
183  const surfaceScalarField Df
184  (
185  dimensionedScalar("viscosity", dimViscosity, 1.0)
186  *diffusivityPtr_->operator()()
187  );
188 
189  volTensorField gradCd("gradCd", fvc::grad(cellDisplacement_));
190 
192 
194  (
196  (
197  2*Df,
198  cellDisplacement_,
199  "laplacian(diffusivity,cellDisplacement)"
200  )
201 
202  + fvc::div
203  (
204  Df
205  *(
207  (
208  cellDisplacement_.mesh().Sf(),
209  gradCd.T() - gradCd
210  )
211 
212  // Solid-body rotation "lambda" term
213  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
214  )
215  )
216 
217  /*
218  - fvc::laplacian
219  (
220  2*Df,
221  cellDisplacement_,
222  "laplacian(diffusivity,cellDisplacement)"
223  )
224 
225  + fvc::div
226  (
227  Df
228  *(
229  fvc::dotInterpolate
230  (
231  cellDisplacement_.mesh().Sf(),
232  gradCd + gradCd.T()
233  )
234 
235  // Solid-body rotation "lambda" term
236  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
237  )
238  )
239  */
240  ==
241  fvOptions(cellDisplacement_)
242  );
243 
246  fvOptions.correct(cellDisplacement_);
247 }
248 
249 
251 (
252  const mapPolyMesh& mpm
253 )
254 {
256 
257  // Update diffusivity. Note two stage to make sure old one is de-registered
258  // before creating/registering new one.
259  diffusivityPtr_.reset(nullptr);
260  diffusivityPtr_ = motionDiffusivity::New
261  (
262  fvMesh_,
263  coeffDict().lookup("diffusivity")
264  );
265 }
266 
267 
268 // ************************************************************************* //
virtual void updateMesh(const mapPolyMesh &)
Update topology.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Virtual base class for displacement motion solver.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
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.
Calculate the matrix for the laplacian of the field.
Finite-volume options.
Definition: fvOptions.H:51
Lookup type of boundary radiation properties.
Definition: lookup.H:57
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:55
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Base class for fvMesh based motionSolvers.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
Calculate the gradient of the given field.
Calculate the laplacian of the given field.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Calculate the divergence of the given field.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Abstract base class for cell-centre mesh motion diffusivity.
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers, to the points. This base class implements the default method which applies volPointInterpolation only.
defineTypeNameAndDebug(combustionModel, 0)
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.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
bool found
A primitive field of type <T> with automated input and output.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127