solidBodyDisplacementLaplacianFvMotionSolver.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) 2016-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "motionInterpolation.H"
30 #include "motionDiffusivity.H"
31 #include "fvmLaplacian.H"
33 #include "OFstream.H"
34 #include "meshTools.H"
35 #include "mapPolyMesh.H"
37 #include "transformField.H"
38 #include "fvOptions.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(solidBodyDisplacementLaplacianFvMotionSolver, 0);
45 
47  (
48  motionSolver,
49  solidBodyDisplacementLaplacianFvMotionSolver,
50  dictionary
51  );
52 
54  (
55  displacementMotionSolver,
56  solidBodyDisplacementLaplacianFvMotionSolver,
57  displacement
58  );
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
64 Foam::solidBodyDisplacementLaplacianFvMotionSolver::
65 solidBodyDisplacementLaplacianFvMotionSolver
66 (
67  const polyMesh& mesh,
68  const IOdictionary& dict
69 )
70 :
73  SBMFPtr_(solidBodyMotionFunction::New(coeffDict(), mesh.time())),
74  cellDisplacement_
75  (
76  IOobject
77  (
78  "cellDisplacement",
79  mesh.time().timeName(),
80  mesh,
81  IOobject::READ_IF_PRESENT,
82  IOobject::AUTO_WRITE
83  ),
84  fvMesh_,
85  dimensionedVector(pointDisplacement_.dimensions(), Zero),
86  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
87  ),
88  pointLocation_(nullptr),
89  interpolationPtr_
90  (
91  coeffDict().found("interpolation")
92  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
93  : motionInterpolation::New(fvMesh_)
94  ),
95  diffusivityPtr_
96  (
97  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
98  ),
99  frozenPointsZone_
100  (
101  coeffDict().found("frozenPointsZone")
102  ? fvMesh_.pointZones().findZoneID
103  (
104  coeffDict().get<word>("frozenPointsZone")
105  )
106  : -1
107  )
108 {
109  IOobject io
110  (
111  "pointLocation",
112  fvMesh_.time().timeName(),
113  fvMesh_,
116  );
117 
118  if (debug)
119  {
120  Info<< "solidBodyDisplacementLaplacianFvMotionSolver:" << nl
121  << " diffusivity : " << diffusivityPtr_().type() << nl
122  << " frozenPoints zone : " << frozenPointsZone_ << endl;
123  }
124 
125 
126  if (io.typeHeaderOk<pointVectorField>(true))
127  {
128  pointLocation_.reset
129  (
130  new pointVectorField
131  (
132  io,
134  )
135  );
136 
137  if (debug)
138  {
139  Info<< "solidBodyDisplacementLaplacianFvMotionSolver :"
140  << " Read pointVectorField "
141  << io.name()
142  << " to be used for boundary conditions on points."
143  << nl
144  << "Boundary conditions:"
145  << pointLocation_().boundaryField().types() << endl;
146  }
147  }
148 }
149 
150 
151 Foam::solidBodyDisplacementLaplacianFvMotionSolver::
152 solidBodyDisplacementLaplacianFvMotionSolver
153 (
154  const polyMesh& mesh,
155  const IOdictionary& dict,
156  const pointVectorField& pointDisplacement,
157  const pointIOField& points0
158 )
159 :
160  displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
161  fvMotionSolver(mesh),
162  SBMFPtr_(solidBodyMotionFunction::New(coeffDict(), mesh.time())),
163  cellDisplacement_
164  (
165  IOobject
166  (
167  "cellDisplacement",
168  mesh.time().timeName(),
169  mesh,
170  IOobject::READ_IF_PRESENT,
171  IOobject::AUTO_WRITE
172  ),
173  fvMesh_,
174  dimensionedVector(pointDisplacement_.dimensions(), Zero),
175  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
176  ),
177  pointLocation_(nullptr),
178  interpolationPtr_
179  (
180  coeffDict().found("interpolation")
181  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
182  : motionInterpolation::New(fvMesh_)
183  ),
184  diffusivityPtr_
185  (
186  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
187  ),
188  frozenPointsZone_
189  (
190  coeffDict().found("frozenPointsZone")
191  ? fvMesh_.pointZones().findZoneID
192  (
193  coeffDict().get<word>("frozenPointsZone")
194  )
195  : -1
196  )
197 {
198  IOobject io
199  (
200  "pointLocation",
201  fvMesh_.time().timeName(),
202  fvMesh_,
205  );
206 
207  if (debug)
208  {
209  Info<< "solidBodyDisplacementLaplacianFvMotionSolver:" << nl
210  << " diffusivity : " << diffusivityPtr_().type() << nl
211  << " frozenPoints zone : " << frozenPointsZone_ << endl;
212  }
213 
214 
215  if (io.typeHeaderOk<pointVectorField>(true))
216  {
217  pointLocation_.reset
218  (
219  new pointVectorField
220  (
221  io,
223  )
224  );
225 
226  if (debug)
227  {
228  Info<< "solidBodyDisplacementLaplacianFvMotionSolver :"
229  << " Read pointVectorField "
230  << io.name()
231  << " to be used for boundary conditions on points."
232  << nl
233  << "Boundary conditions:"
234  << pointLocation_().boundaryField().types() << endl;
235  }
236  }
237 }
238 
239 
240 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
241 
244 {}
245 
246 
247 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
248 
251 {
252  if (!diffusivityPtr_)
253  {
254  diffusivityPtr_ = motionDiffusivity::New
255  (
256  fvMesh_,
257  coeffDict().lookup("diffusivity")
258  );
259  }
260 
261  return *diffusivityPtr_;
262 }
263 
264 
267 {
268  interpolationPtr_->interpolate
269  (
270  cellDisplacement_,
271  pointDisplacement_
272  );
273 
274  tmp<pointField> tnewPoints
275  (
276  transformPoints(SBMFPtr_().transformation(), points0())
277  );
278  const pointField& newPoints = tnewPoints();
279 
280  if (pointLocation_)
281  {
282  if (debug)
283  {
284  Info<< "solidBodyDisplacementLaplacianFvMotionSolver : applying "
285  << " boundary conditions on " << pointLocation_().name()
286  << " to new point location."
287  << endl;
288  }
289 
290  pointLocation_().primitiveFieldRef() =
291  newPoints
292  + pointDisplacement_.internalField();
293 
294  pointLocation_().correctBoundaryConditions();
295 
296  // Implement frozen points
297  if (frozenPointsZone_ != -1)
298  {
299  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
300 
301  forAll(pz, i)
302  {
303  pointLocation_()[pz[i]] = newPoints[pz[i]];
304  }
305  }
306 
307  twoDCorrectPoints(pointLocation_().primitiveFieldRef());
308 
309  return tmp<pointField>(pointLocation_().primitiveField());
310  }
311  else
312  {
313  tmp<pointField> tcurPoints
314  (
315  newPoints + pointDisplacement_.primitiveField()
316  );
317  pointField& curPoints = tcurPoints.ref();
318 
319  // Implement frozen points
320  if (frozenPointsZone_ != -1)
321  {
322  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
323 
324  forAll(pz, i)
325  {
326  curPoints[pz[i]] = newPoints[pz[i]];
327  }
328  }
329 
330  twoDCorrectPoints(curPoints);
331 
332  return tcurPoints;
333  }
334 }
335 
336 
338 {
339  // The points have moved so before interpolation update
340  // the motionSolver accordingly
341  movePoints(fvMesh_.points());
342 
343  diffusivity().correct();
344  pointDisplacement_.boundaryFieldRef().updateCoeffs();
345 
346  fv::options& fvOptions(fv::options::New(fvMesh_));
347 
349  (
351  (
352  dimensionedScalar("viscosity", dimViscosity, 1.0)
353  *diffusivity().operator()(),
354  cellDisplacement_,
355  "laplacian(diffusivity,cellDisplacement)"
356  )
357  ==
358  fvOptions(cellDisplacement_)
359  );
360 
363  fvOptions.correct(cellDisplacement_);
364 }
365 
366 
368 (
369  const mapPolyMesh& mpm
370 )
371 {
373 
374  // Update diffusivity. Note two stage to make sure old one is de-registered
375  // before creating/registering new one.
376  diffusivityPtr_.clear();
377 }
378 
379 
380 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
dictionary dict
Virtual base class for displacement motion solver.
void transformPoints(vectorField &, const septernion &, const vectorField &)
Transform given vectorField of coordinates with the given septernion.
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:38
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const dimensionSet dimViscosity
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
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.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
const fvMesh & fvMesh_
The fvMesh to be moved.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Lookup type of boundary radiation properties.
Definition: lookup.H:57
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:55
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
conserve primitiveFieldRef()+
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Base class for fvMesh based motionSolvers.
Spatial transformation functions for primitive fields.
word timeName
Definition: getTimeIndex.H:3
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
Base class for defining solid-body motions.
motionDiffusivity & diffusivity()
Return reference to the diffusivity field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Vector< scalar > vector
Definition: vector.H:57
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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)
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
IOobject(const IOobject &)=default
Copy construct.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
Automatically write from objectRegistry::writeObject()
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
bool found
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127