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-2024 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 :
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  diffusivityType_(dict.lookup("diffusivity")),
83  diffusivityPtr_
84  (
85  motionDiffusivity::New(fvMesh_, diffusivityType_)
86  )
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91 
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
101 {
103  (
104  cellDisplacement_,
105  pointDisplacement_
106  );
107 
108  tmp<pointField> tcurPoints
109  (
110  points0() + pointDisplacement().primitiveField()
111  );
112 
113  twoDCorrectPoints(tcurPoints.ref());
114 
115  return tcurPoints;
116 }
117 
118 
120 {
121  // The points have moved so before interpolation update
122  // the mtionSolver accordingly
123  movePoints(fvMesh_.points());
124 
125  diffusivityPtr_->correct();
126  pointDisplacement_.boundaryFieldRef().updateCoeffs();
127 
128  surfaceScalarField Df(diffusivityPtr_->operator()());
129 
130  volTensorField gradCd("gradCd", fvc::grad(cellDisplacement_));
131 
133  (
135  (
136  2*Df,
137  cellDisplacement_,
138  "laplacian(diffusivity,cellDisplacement)"
139  )
140 
141  + fvc::div
142  (
143  Df
144  *(
146  (
147  cellDisplacement_.mesh().Sf(),
148  gradCd.T() - gradCd
149  )
150 
151  // Solid-body rotation "lambda" term
152  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
153  )
154  )
155 
156  /*
157  - fvc::laplacian
158  (
159  2*Df,
160  cellDisplacement_,
161  "laplacian(diffusivity,cellDisplacement)"
162  )
163 
164  + fvc::div
165  (
166  Df
167  *(
168  fvc::dotInterpolate
169  (
170  cellDisplacement_.mesh().Sf(),
171  gradCd + gradCd.T()
172  )
173 
174  // Solid-body rotation "lambda" term
175  - cellDisplacement_.mesh().Sf()*fvc::interpolate(tr(gradCd))
176  )
177  )
178  */
179  );
180 }
181 
182 
184 (
185  const polyTopoChangeMap& map
186 )
187 {
189 
190  // Update diffusivity. Note two stage to make sure old one is de-registered
191  // before creating/registering new one.
192  diffusivityPtr_.reset(nullptr);
193  diffusivityType_.rewind();
194  diffusivityPtr_ = motionDiffusivity::New
195  (
196  fvMesh_,
197  diffusivityType_
198  );
199 }
200 
201 
203 (
204  const polyMeshMap& map
205 )
206 {
208 
209  // Update diffusivity. Note two stage to make sure old one is de-registered
210  // before creating/registering new one.
211  diffusivityPtr_.reset(nullptr);
212  diffusivityType_.rewind();
213  diffusivityPtr_ = motionDiffusivity::New
214  (
215  fvMesh_,
216  diffusivityType_
217  );
218 }
219 
220 
221 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static volPointInterpolation & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricField class.
tmp< GeometricField< Type, GeoMesh, Field > > 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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:197
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
defineTypeNameAndDebug(combustionModel, 0)
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dictionary dict