elasticityMotionSolver.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) 2007-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "elasticityMotionSolver.H"
31 #include "motionInterpolation.H"
32 #include "wallDist.H"
34 #include "fvMatrices.H"
35 #include "fvcDiv.H"
36 #include "fvmDiv.H"
37 #include "fvmDiv.H"
38 #include "fvmLaplacian.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(elasticityMotionSolver, 1);
46 
48  (
49  motionSolver,
50  elasticityMotionSolver,
51  dictionary
52  );
53 }
54 
55 
56 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
57 
59 {
60  // Adjust boundary conditions based on the steps to be executed
62  {
63  pointPatchVectorField& pointBCs =
65  if (isA<fixedValuePointPatchVectorField>(pointBCs))
66  {
67  auto& fixedValueBCs =
68  refCast<fixedValuePointPatchVectorField>(pointBCs);
69  fixedValueBCs == fixedValueBCs/scalar(nSteps_);
70  }
71  }
72 
73  // Copy boundary conditions to internalField
74  // Needed for interpolation to faces
75  pointMotionU_.boundaryFieldRef().updateCoeffs();
76 
77  // Interpolate boundary conditions from points to faces
79  {
81  if (isA<fixedValueFvPatchVectorField>(bField))
82  {
83  const pointField& points = fvMesh_.points();
84  const polyPatch& patch = mesh().boundaryMesh()[pI];
85  forAll(bField, fI)
86  {
87  bField[fI] = patch[fI].average(points, pointMotionU_);
88  }
89  }
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
96 Foam::elasticityMotionSolver::elasticityMotionSolver
97 (
98  const polyMesh& mesh,
99  const IOdictionary& dict
100 )
101 :
102  motionSolver(mesh, dict, typeName),
103  fvMesh_
104  (
105  const_cast<fvMesh&>
106  (
107  refCast<const fvMesh>(mesh)
108  )
109  ),
110  pointMotionU_
111  (
112  IOobject
113  (
114  "pointMotionU",
115  mesh.time().timeName(),
116  mesh,
119  ),
122  fixedValuePointPatchVectorField::typeName
123  ),
124  cellMotionU_
125  (
126  IOobject
127  (
128  "cellMotionU",
129  mesh.time().timeName(),
130  mesh,
133  ),
134  fvMesh_,
135  dimensionedVector(pointMotionU_.dimensions(), Zero),
136  pointMotionU_.boundaryField().types()
137  ),
138  interpolationPtr_
139  (
140  coeffDict().found("interpolation")
141  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
142  : motionInterpolation::New(fvMesh_)
143  ),
144  E_
145  (
146  IOobject
147  (
148  "mu",
149  mesh.time().timeName(),
150  mesh,
151  IOobject::NO_READ,
152  IOobject::NO_WRITE
153  ),
154  fvMesh_,
156  zeroGradientFvPatchScalarField::typeName
157  ),
158  exponent_(this->coeffDict().get<scalar>("exponent")),
159  nSteps_(this->coeffDict().get<label>("steps")),
160  nIters_(this->coeffDict().get<label>("iters")),
161  tolerance_(this->coeffDict().get<scalar>("tolerance"))
162 {}
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
168 {
169  tmp<pointField> tnewPoints(new pointField(mesh().points()));
170 
171  return tnewPoints;
172 }
173 
174 
176 {
177  // Re-init to zero
178  cellMotionU_.primitiveFieldRef() = vector::zero;
179 
180  // Adjust boundary conditions based on the number of steps to be executed
181  // and interpolate to faces
182  setBoundaryConditions();
183 
184  // Solve the elasticity equations in a stepped manner
185  for (label istep = 0; istep < nSteps_; ++istep)
186  {
187  Info<< "Step " << istep << endl;
188 
189  // Update diffusivity
190  const scalarField& vols = mesh().cellVolumes();
191  E_.primitiveFieldRef() = 1./pow(vols, exponent_);
192  E_.correctBoundaryConditions();
193 
194  for (label iter = 0; iter < nIters_; ++iter)
195  {
196  Info<< "Iteration " << iter << endl;
197  cellMotionU_.storePrevIter();
198  fvVectorMatrix dEqn
199  (
200  fvm::laplacian(2*E_, cellMotionU_)
201  + fvc::div(2*E_*T(fvc::grad(cellMotionU_)))
202  - fvc::div(E_*fvc::div(cellMotionU_)*tensor::I)
203  );
204 
205  scalar residual = mag(dEqn.solve().initialResidual());
206  cellMotionU_.relax();
207 
208  // Print execution time
209  fvMesh_.time().printExecutionTime(Info);
210 
211  // Check convergence
212  if (residual < tolerance_)
213  {
214  Info<< "\n***Reached mesh movement convergence limit for step "
215  << istep
216  << " iteration " << iter << "***\n\n";
217  break;
218  }
219  }
220 
221  // Interpolate from cells to points
222  interpolationPtr_->interpolate(cellMotionU_, pointMotionU_);
223  vectorField newPoints
224  (
225  mesh().points() + pointMotionU_.primitiveFieldRef()
226  );
227 
228  // Move points and check mesh
229  fvMesh_.movePoints(newPoints);
230  fvMesh_.checkMesh(true);
231 
232  if (debug)
233  {
234  Info<< " Writing new mesh points " << endl;
236  (
237  IOobject
238  (
239  "points",
240  mesh().pointsInstance(),
241  mesh().meshSubDir,
242  mesh(),
245  ),
246  mesh().points()
247  );
248  points.write();
249  }
250  }
251 }
252 
253 
255 {
256  // Do nothing
257 }
258 
259 
261 {
262  // Do nothing
263 }
264 
265 
266 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::fvPatchField< vector >
wallDist.H
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::Tensor::I
static const Tensor I
Definition: Tensor.H:85
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::elasticityMotionSolver::solve
virtual void solve()
Solve for motion.
Definition: elasticityMotionSolver.C:175
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
fvcDiv.H
Calculate the divergence of the given field.
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::elasticityMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location. Mesh is actually moved in solve()
Definition: elasticityMotionSolver.C:167
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::pointPatchField< vector >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::elasticityMotionSolver::cellMotionU_
volVectorField cellMotionU_
Definition: elasticityMotionSolver.H:81
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::elasticityMotionSolver::fvMesh_
fvMesh & fvMesh_
Definition: elasticityMotionSolver.H:79
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:96
fixedValuePointPatchFields.H
Foam::elasticityMotionSolver::movePoints
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Definition: elasticityMotionSolver.C:254
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::motionInterpolation
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
Definition: motionInterpolation.H:52
Foam::motionSolver
Virtual base class for mesh motion solver.
Definition: motionSolver.H:58
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::motionSolver::mesh
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:144
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::elasticityMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update the mesh corresponding to given map.
Definition: elasticityMotionSolver.C:260
Foam::elasticityMotionSolver::setBoundaryConditions
void setBoundaryConditions()
Set boundary conditions of cellMotionU based on pointMotionU.
Definition: elasticityMotionSolver.C:58
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::motionInterpolation::New
static autoPtr< motionInterpolation > New(const fvMesh &mesh)
Select default.
Definition: motionInterpolation.C:66
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
motionInterpolation.H
Foam::elasticityMotionSolver::pointMotionU_
pointVectorField pointMotionU_
Definition: elasticityMotionSolver.H:80
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
elasticityMotionSolver.H
Foam::elasticityMotionSolver::nSteps_
label nSteps_
Intermediate steps to solve the PDEs.
Definition: elasticityMotionSolver.H:93