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 | Copyright (C) 2011-2017 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "motionDiffusivity.H"
28 #include "fvmLaplacian.H"
30 #include "fvcDiv.H"
31 #include "fvcGrad.H"
32 #include "surfaceInterpolate.H"
33 #include "fvcLaplacian.H"
34 #include "mapPolyMesh.H"
35 #include "volPointInterpolation.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(displacementSBRStressFvMotionSolver, 0);
42 
44  (
45  motionSolver,
46  displacementSBRStressFvMotionSolver,
47  dictionary
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 Foam::displacementSBRStressFvMotionSolver::displacementSBRStressFvMotionSolver
55 (
56  const polyMesh& mesh,
57  const IOdictionary& dict
58 )
59 :
60  displacementMotionSolver(mesh, dict, typeName),
61  fvMotionSolver(mesh),
62  cellDisplacement_
63  (
64  IOobject
65  (
66  "cellDisplacement",
67  mesh.time().timeName(),
68  mesh,
71  ),
72  fvMesh_,
74  (
75  "cellDisplacement",
76  pointDisplacement().dimensions(),
77  Zero
78  ),
79  cellMotionBoundaryTypes<vector>(pointDisplacement().boundaryField())
80  ),
81  diffusivityPtr_
82  (
83  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
84  )
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
99 {
101  (
102  cellDisplacement_,
103  pointDisplacement_
104  );
105 
106  tmp<pointField> tcurPoints
107  (
108  points0() + pointDisplacement().primitiveField()
109  );
110 
111  twoDCorrectPoints(tcurPoints.ref());
112 
113  return tcurPoints;
114 }
115 
116 
118 {
119  // The points have moved so before interpolation update
120  // the mtionSolver accordingly
121  movePoints(fvMesh_.points());
122 
123  diffusivityPtr_->correct();
124  pointDisplacement_.boundaryFieldRef().updateCoeffs();
125 
126  surfaceScalarField Df(diffusivityPtr_->operator()());
127 
128  volTensorField gradCd("gradCd", fvc::grad(cellDisplacement_));
129 
131  (
133  (
134  2*Df,
135  cellDisplacement_,
136  "laplacian(diffusivity,cellDisplacement)"
137  )
138 
139  + fvc::div
140  (
141  Df
142  *(
144  (
145  cellDisplacement_.mesh().Sf(),
146  gradCd.T() - gradCd
147  )
148 
149  // Solid-body rotation "lambda" term
150  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
151  )
152  )
153 
154  /*
155  - fvc::laplacian
156  (
157  2*Df,
158  cellDisplacement_,
159  "laplacian(diffusivity,cellDisplacement)"
160  )
161 
162  + fvc::div
163  (
164  Df
165  *(
166  fvc::dotInterpolate
167  (
168  cellDisplacement_.mesh().Sf(),
169  gradCd + gradCd.T()
170  )
171 
172  // Solid-body rotation "lambda" term
173  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
174  )
175  )
176  */
177  );
178 }
179 
180 
182 (
183  const mapPolyMesh& mpm
184 )
185 {
187 
188  // Update diffusivity. Note two stage to make sure old one is de-registered
189  // before creating/registering new one.
190  diffusivityPtr_.reset(nullptr);
191  diffusivityPtr_ = motionDiffusivity::New
192  (
193  fvMesh_,
194  coeffDict().lookup("diffusivity")
195  );
196 }
197 
198 
199 // ************************************************************************* //
virtual void updateMesh(const mapPolyMesh &)
Update topology.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Virtual base class for displacement motion solver.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.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.
Calculate the matrix for the laplacian of the field.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
static const volPointInterpolation & New(const fvMesh &mesh)
Definition: MeshObject.C:44
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
Base class for fvMesh based motionSolvers.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Calculate the gradient of the given field.
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Calculate the laplacian of the given field.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
static const zero Zero
Definition: zero.H:91
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Calculate the divergence of the given field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Time & time() const
Return time.
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.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
A class for managing temporary objects.
Definition: PtrList.H:53
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.