displacementSBRStressFvMotionSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 "polyTopoChangeMap.H"
35 #include "volPointInterpolation.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
42 
44  (
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const word& name,
57  const polyMesh& mesh,
58  const dictionary& dict
59 )
60 :
61  displacementMotionSolver(name, mesh, dict, typeName),
62  fvMotionSolver(mesh),
63  cellDisplacement_
64  (
65  IOobject
66  (
67  "cellDisplacement",
68  mesh.time().name(),
69  mesh,
70  IOobject::READ_IF_PRESENT,
71  IOobject::AUTO_WRITE
72  ),
73  fvMesh_,
75  (
76  "cellDisplacement",
77  pointDisplacement().dimensions(),
78  Zero
79  ),
80  cellMotionBoundaryTypes<vector>(pointDisplacement().boundaryField())
81  ),
82  diffusivityPtr_
83  (
84  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
85  )
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
100 {
102  (
103  cellDisplacement_,
104  pointDisplacement_
105  );
106 
107  tmp<pointField> tcurPoints
108  (
109  points0() + pointDisplacement().primitiveField()
110  );
111 
112  twoDCorrectPoints(tcurPoints.ref());
113 
114  return tcurPoints;
115 }
116 
117 
119 {
120  // The points have moved so before interpolation update
121  // the mtionSolver accordingly
122  movePoints(fvMesh_.points());
123 
124  diffusivityPtr_->correct();
125  pointDisplacement_.boundaryFieldRef().updateCoeffs();
126 
127  surfaceScalarField Df(diffusivityPtr_->operator()());
128 
129  volTensorField gradCd("gradCd", fvc::grad(cellDisplacement_));
130 
132  (
134  (
135  2*Df,
136  cellDisplacement_,
137  "laplacian(diffusivity,cellDisplacement)"
138  )
139 
140  + fvc::div
141  (
142  Df
143  *(
145  (
146  cellDisplacement_.mesh().Sf(),
147  gradCd.T() - gradCd
148  )
149 
150  // Solid-body rotation "lambda" term
151  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
152  )
153  )
154 
155  /*
156  - fvc::laplacian
157  (
158  2*Df,
159  cellDisplacement_,
160  "laplacian(diffusivity,cellDisplacement)"
161  )
162 
163  + fvc::div
164  (
165  Df
166  *(
167  fvc::dotInterpolate
168  (
169  cellDisplacement_.mesh().Sf(),
170  gradCd + gradCd.T()
171  )
172 
173  // Solid-body rotation "lambda" term
174  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
175  )
176  )
177  */
178  );
179 }
180 
181 
183 (
184  const polyTopoChangeMap& map
185 )
186 {
188 
189  // Update diffusivity. Note two stage to make sure old one is de-registered
190  // before creating/registering new one.
191  diffusivityPtr_.reset(nullptr);
192  diffusivityPtr_ = motionDiffusivity::New
193  (
194  fvMesh_,
195  coeffDict().lookup("diffusivity")
196  );
197 }
198 
199 
201 (
202  const polyMeshMap& map
203 )
204 {
206 
207  // Update diffusivity. Note two stage to make sure old one is de-registered
208  // before creating/registering new one.
209  diffusivityPtr_.reset(nullptr);
210  diffusivityPtr_ = motionDiffusivity::New
211  (
212  fvMesh_,
213  coeffDict().lookup("diffusivity")
214  );
215 }
216 
217 
218 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Virtual base class for displacement motion solver.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Mesh motion solver for an fvMesh. Based on solving the cell-centre solid-body rotation stress equatio...
displacementSBRStressFvMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Base class for fvMesh based motionSolvers.
Abstract base class for cell-centre mesh motion diffusivity.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
A class for handling words, derived from string.
Definition: word.H:62
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.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dictionary dict