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-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 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 "motionDiffusivity.H"
33 #include "wallDist.H"
35 #include "fvMatrices.H"
36 #include "fvcDiv.H"
37 #include "fvmDiv.H"
38 #include "fvmDiv.H"
39 #include "fvmLaplacian.H"
40 #include "surfaceInterpolate.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  defineTypeNameAndDebug(elasticityMotionSolver, 1);
48 
50  (
51  motionSolver,
52  elasticityMotionSolver,
54  );
55 }
56 
57 
58 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
59 
61 {
62  // Adjust boundary conditions based on the steps to be executed
64  {
65  fvPatchVectorField& bc =
67  if (isA<fixedValueFvPatchVectorField>(bc))
68  {
69  auto& fixedValueBCs =
70  refCast<fixedValueFvPatchVectorField>(bc);
71  fixedValueBCs == fixedValueBCs/scalar(nSteps_);
72  }
73  }
74  /*
75  // Adjust boundary conditions based on the steps to be executed
76  forAll(pointMotionU_.boundaryField(), patchI)
77  {
78  pointPatchVectorField& pointBCs =
79  pointMotionU_.boundaryFieldRef()[patchI];
80  if (isA<fixedValuePointPatchVectorField>(pointBCs))
81  {
82  auto& fixedValueBCs =
83  refCast<fixedValuePointPatchVectorField>(pointBCs);
84  fixedValueBCs == fixedValueBCs/scalar(nSteps_);
85  }
86  }
87 
88  // Copy boundary conditions to internalField
89  // Needed for interpolation to faces
90  pointMotionU_.boundaryFieldRef().updateCoeffs();
91 
92  // Interpolate boundary conditions from points to faces
93  forAll(cellMotionU_.boundaryField(), pI)
94  {
95  fvPatchVectorField& bField = cellMotionU_.boundaryFieldRef()[pI];
96  if (isA<fixedValueFvPatchVectorField>(bField))
97  {
98  const pointField& points = fvMesh_.points();
99  const polyPatch& patch = mesh().boundaryMesh()[pI];
100  forAll(bField, fI)
101  {
102  bField[fI] = patch[fI].average(points, pointMotionU_);
103  }
104  }
105  }
106  */
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111 
112 Foam::elasticityMotionSolver::elasticityMotionSolver
113 (
114  const polyMesh& mesh,
115  const IOdictionary& dict
116 )
117 :
118  motionSolver(mesh, dict, typeName),
119  fvMesh_
120  (
121  const_cast<fvMesh&>
122  (
123  refCast<const fvMesh>(mesh)
124  )
125  ),
126  pointMotionU_
127  (
128  IOobject
129  (
130  "pointMotionU",
131  mesh.time().timeName(),
132  mesh,
133  IOobject::READ_IF_PRESENT,
134  IOobject::AUTO_WRITE
135  ),
136  pointMesh::New(mesh),
138  fixedValuePointPatchVectorField::typeName
139  ),
140  cellMotionU_
141  (
142  IOobject
143  (
144  "cellMotionU",
145  mesh.time().timeName(),
146  mesh,
147  IOobject::READ_IF_PRESENT,
148  IOobject::AUTO_WRITE
149  ),
150  fvMesh_,
151  dimensionedVector(pointMotionU_.dimensions(), Zero),
152  pointMotionU_.boundaryField().types()
153  ),
154  interpolationPtr_
155  (
156  coeffDict().found("interpolation")
157  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
158  : motionInterpolation::New(fvMesh_)
159  ),
160  diffusivityPtr_
161  (
162  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
163  ),
164  nSteps_(this->coeffDict().get<label>("steps")),
165  nIters_(this->coeffDict().get<label>("iters")),
166  tolerance_(this->coeffDict().get<scalar>("tolerance"))
167 {}
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
173 {
174  tmp<pointField> tnewPoints(new pointField(mesh().points()));
175 
176  return tnewPoints;
177 }
178 
179 
181 {
182  // Re-init to zero
183  cellMotionU_.primitiveFieldRef() = Zero;
184 
185  // Adjust boundary conditions based on the number of steps to be executed
186  // and interpolate to faces
187  setBoundaryConditions();
188 
189  // Solve the elasticity equations in a stepped manner
190  for (label istep = 0; istep < nSteps_; ++istep)
191  {
192  Info<< "Step " << istep << endl;
193 
194  // Update diffusivity
195  diffusivityPtr_->correct();
196  const surfaceScalarField E(diffusivityPtr_->operator()());
197  const surfaceVectorField& Sf = fvMesh_.Sf();
198  for (label iter = 0; iter < nIters_; ++iter)
199  {
200  Info<< "Iteration " << iter << endl;
201  cellMotionU_.storePrevIter();
202  fvVectorMatrix dEqn
203  (
204  fvm::laplacian(2*E, cellMotionU_)
205  + fvc::div(2*E*(fvc::interpolate(fvc::grad(cellMotionU_)) & Sf))
206  - fvc::div(E*fvc::interpolate(fvc::div(cellMotionU_))*Sf)
207  );
208 
209  scalar residual = mag(dEqn.solve().initialResidual());
210  cellMotionU_.relax();
211 
212  // Print execution time
213  fvMesh_.time().printExecutionTime(Info);
214 
215  // Check convergence
216  if (residual < tolerance_)
217  {
218  Info<< "\n***Reached mesh movement convergence limit for step "
219  << istep
220  << " iteration " << iter << "***\n\n";
221  break;
222  }
223  }
224 
225  interpolationPtr_->interpolate
226  (
227  cellMotionU_,
228  pointMotionU_
229  );
230 
231  tmp<vectorField> newPoints
232  (
233  fvMesh_.points() + pointMotionU_.primitiveField()
234  );
235 
236  /*
237  // Interpolate from cells to points
238  interpolationPtr_->interpolate(cellMotionU_, pointMotionU_);
239 
240  syncTools::syncPointList
241  (
242  fvMesh_,
243  pointMotionU_.primitiveFieldRef(),
244  maxEqOp<vector>(),
245  vector::zero
246  );
247 
248  vectorField newPoints
249  (
250  mesh().points() + pointMotionU_.primitiveFieldRef()
251  );
252  */
253 
254  // Move points and check mesh
255  fvMesh_.movePoints(newPoints);
256  fvMesh_.checkMesh(true);
257 
258  if (debug)
259  {
260  Info<< " Writing new mesh points " << endl;
262  (
263  IOobject
264  (
265  "points",
266  mesh().pointsInstance(),
267  mesh().meshSubDir,
268  mesh(),
271  ),
272  mesh().points()
273  );
274  points.write();
275  }
276  }
277 }
278 
281 {
282  // Do nothing
283 }
284 
285 
287 {
288  // Update diffusivity. Note two stage to make sure old one is de-registered
289  // before creating/registering new one.
290  diffusivityPtr_.reset(nullptr);
291  diffusivityPtr_ = motionDiffusivity::New
292  (
293  fvMesh_,
294  coeffDict().lookup("diffusivity")
295  );
296 }
297 
298 
299 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
label nSteps_
Intermediate steps to solve the PDEs.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:38
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
virtual void solve()
Solve for motion.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:54
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
void setBoundaryConditions()
Set boundary conditions of cellMotionU based on pointMotionU.
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.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
const pointField & points
virtual void updateMesh(const mapPolyMesh &)
Update the mesh corresponding to given map.
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
Calculate the divergence of the given field.
Abstract base class for cell-centre mesh motion diffusivity.
int debug
Static debugging option.
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.
Calculate the matrix for the divergence of the given field and flux.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
bool found
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
virtual tmp< pointField > curPoints() const
Return point location. Mesh is actually moved in solve()
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127