displacementLaplacianFvMotionSolver.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "motionInterpolation.H"
31 #include "motionDiffusivity.H"
32 #include "fvmLaplacian.H"
34 #include "OFstream.H"
35 #include "meshTools.H"
36 #include "mapPolyMesh.H"
37 #include "fvOptions.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(displacementLaplacianFvMotionSolver, 0);
44 
46  (
47  motionSolver,
48  displacementLaplacianFvMotionSolver,
49  dictionary
50  );
51 
53  (
54  displacementMotionSolver,
55  displacementLaplacianFvMotionSolver,
56  displacement
57  );
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
63 Foam::displacementLaplacianFvMotionSolver::displacementLaplacianFvMotionSolver
64 (
65  const polyMesh& mesh,
66  const IOdictionary& dict
67 )
68 :
71  cellDisplacement_
72  (
73  IOobject
74  (
75  "cellDisplacement",
76  mesh.time().timeName(),
77  mesh,
78  IOobject::READ_IF_PRESENT,
79  IOobject::AUTO_WRITE
80  ),
81  fvMesh_,
82  dimensionedVector(pointDisplacement_.dimensions(), Zero),
83  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
84  ),
85  pointLocation_(nullptr),
86  interpolationPtr_
87  (
88  coeffDict().found("interpolation")
89  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
90  : motionInterpolation::New(fvMesh_)
91  ),
92  diffusivityPtr_
93  (
94  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
95  ),
96  frozenPointsZone_
97  (
98  coeffDict().found("frozenPointsZone")
99  ? fvMesh_.pointZones().findZoneID
100  (
101  coeffDict().get<word>("frozenPointsZone")
102  )
103  : -1
104  )
105 {
106  IOobject io
107  (
108  "pointLocation",
109  fvMesh_.time().timeName(),
110  fvMesh_,
111  IOobject::MUST_READ,
112  IOobject::AUTO_WRITE
113  );
114 
115  if (debug)
116  {
117  Info<< "displacementLaplacianFvMotionSolver:" << nl
118  << " diffusivity : " << diffusivityPtr_().type() << nl
119  << " frozenPoints zone : " << frozenPointsZone_ << endl;
120  }
121 
122 
123  if (io.typeHeaderOk<pointVectorField>(true))
124  {
125  pointLocation_.reset
126  (
127  new pointVectorField
128  (
129  io,
130  pointMesh::New(fvMesh_)
131  )
132  );
133 
134  if (debug)
135  {
136  Info<< "displacementLaplacianFvMotionSolver :"
137  << " Read pointVectorField "
138  << io.name()
139  << " to be used for boundary conditions on points."
140  << nl
141  << "Boundary conditions:"
142  << pointLocation_().boundaryField().types() << endl;
143  }
144  }
145 }
146 
147 
148 Foam::displacementLaplacianFvMotionSolver::
149 displacementLaplacianFvMotionSolver
150 (
151  const polyMesh& mesh,
152  const IOdictionary& dict,
153  const pointVectorField& pointDisplacement,
154  const pointIOField& points0
155 )
156 :
157  displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
159  cellDisplacement_
160  (
161  IOobject
162  (
163  "cellDisplacement",
164  mesh.time().timeName(),
165  mesh,
166  IOobject::READ_IF_PRESENT,
167  IOobject::AUTO_WRITE
168  ),
169  fvMesh_,
170  dimensionedVector(pointDisplacement_.dimensions(), Zero),
171  cellMotionBoundaryTypes<vector>(pointDisplacement_.boundaryField())
172  ),
173  pointLocation_(nullptr),
174  interpolationPtr_
175  (
176  coeffDict().found("interpolation")
177  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
178  : motionInterpolation::New(fvMesh_)
179  ),
180  diffusivityPtr_
181  (
182  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
183  ),
184  frozenPointsZone_
185  (
186  coeffDict().found("frozenPointsZone")
187  ? fvMesh_.pointZones().findZoneID
188  (
189  coeffDict().get<word>("frozenPointsZone")
190  )
191  : -1
192  )
193 {
194  IOobject io
195  (
196  "pointLocation",
197  fvMesh_.time().timeName(),
198  fvMesh_,
199  IOobject::MUST_READ,
200  IOobject::AUTO_WRITE
201  );
202 
203  if (debug)
204  {
205  Info<< "displacementLaplacianFvMotionSolver:" << nl
206  << " diffusivity : " << diffusivityPtr_().type() << nl
207  << " frozenPoints zone : " << frozenPointsZone_ << endl;
208  }
209 
210 
211  if (io.typeHeaderOk<pointVectorField>(true))
212  {
213  pointLocation_.reset
214  (
215  new pointVectorField
216  (
217  io,
218  pointMesh::New(fvMesh_)
219  )
220  );
221 
222  if (debug)
223  {
224  Info<< "displacementLaplacianFvMotionSolver :"
225  << " Read pointVectorField "
226  << io.name()
227  << " to be used for boundary conditions on points."
228  << nl
229  << "Boundary conditions:"
230  << pointLocation_().boundaryField().types() << endl;
231  }
232  }
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
237 
240 {}
241 
242 
243 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244 
247 {
248  if (!diffusivityPtr_)
249  {
250  diffusivityPtr_ = motionDiffusivity::New
251  (
252  fvMesh_,
253  coeffDict().lookup("diffusivity")
254  );
255  }
256 
257  return *diffusivityPtr_;
258 }
259 
260 
263 {
264  interpolationPtr_->interpolate
265  (
266  cellDisplacement_,
267  pointDisplacement_
268  );
269 
270  if (pointLocation_)
271  {
272  if (debug)
273  {
274  Info<< "displacementLaplacianFvMotionSolver : applying "
275  << " boundary conditions on " << pointLocation_().name()
276  << " to new point location."
277  << endl;
278  }
279 
280  pointLocation_().primitiveFieldRef() =
281  points0()
282  + pointDisplacement_.primitiveField();
283 
284  pointLocation_().correctBoundaryConditions();
285 
286  // Implement frozen points
287  if (frozenPointsZone_ != -1)
288  {
289  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
290 
291  forAll(pz, i)
292  {
293  pointLocation_()[pz[i]] = points0()[pz[i]];
294  }
295  }
296 
297  twoDCorrectPoints(pointLocation_().primitiveFieldRef());
298 
299  return tmp<pointField>(pointLocation_().primitiveField());
300  }
301  else
302  {
303  tmp<pointField> tcurPoints
304  (
305  points0() + pointDisplacement_.primitiveField()
306  );
307  pointField& curPoints = tcurPoints.ref();
308 
309  // Implement frozen points
310  if (frozenPointsZone_ != -1)
311  {
312  const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
313 
314  forAll(pz, i)
315  {
316  curPoints[pz[i]] = points0()[pz[i]];
317  }
318  }
319 
320  twoDCorrectPoints(curPoints);
321 
322  return tcurPoints;
323  }
324 }
325 
326 
328 {
329  // The points have moved so before interpolation update
330  // the motionSolver accordingly
331  movePoints(fvMesh_.points());
332 
333  diffusivity().correct();
334  pointDisplacement_.boundaryFieldRef().updateCoeffs();
335 
337 
338  // We explicitly do NOT want to interpolate the motion inbetween
339  // different regions so bypass all the matrix manipulation.
341  (
343  (
344  dimensionedScalar("viscosity", dimViscosity, 1.0)
345  *diffusivity().operator()(),
346  cellDisplacement_,
347  "laplacian(diffusivity,cellDisplacement)"
348  )
349  ==
350  fvOptions(cellDisplacement_)
351  );
352 
355  fvOptions.correct(cellDisplacement_);
356 }
357 
358 
360 (
361  const mapPolyMesh& mpm
362 )
363 {
365 
366  // Update diffusivity. Note two stage to make sure old one is de-registered
367  // before creating/registering new one.
368  diffusivityPtr_.clear();
369 }
370 
371 
372 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::displacementLaplacianFvMotionSolver::~displacementLaplacianFvMotionSolver
~displacementLaplacianFvMotionSolver()
Destructor.
Definition: displacementLaplacianFvMotionSolver.C:239
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
meshTools.H
fvOptions.H
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::displacementLaplacianFvMotionSolver::solve
virtual void solve()
Solve for motion.
Definition: displacementLaplacianFvMotionSolver.C:327
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::motionDiffusivity
Abstract base class for cell-centre mesh motion diffusivity.
Definition: motionDiffusivity.H:51
Foam::fv::optionList::constrain
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
Definition: fvOptionListTemplates.C:314
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::fvMatrix::solveSegregatedOrCoupled
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:62
mapPolyMesh.H
Foam::pointZone
A subset of mesh points.
Definition: pointZone.H:65
displacementLaplacianFvMotionSolver.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::displacementLaplacianFvMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update topology.
Definition: displacementLaplacianFvMotionSolver.C:360
primitiveFieldRef
conserve primitiveFieldRef()+
Foam::displacementLaplacianFvMotionSolver::diffusivity
motionDiffusivity & diffusivity()
Return reference to the diffusivity field.
Definition: displacementLaplacianFvMotionSolver.C:246
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
TEqn
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
Foam::fvMotionSolver
Base class for fvMesh based motionSolvers.
Definition: fvMotionSolver.H:50
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
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::motionDiffusivity::New
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Definition: motionDiffusivity.C:51
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dimViscosity
const dimensionSet dimViscosity
Foam::displacementMotionSolver
Virtual base class for displacement motion solver.
Definition: displacementMotionSolver.H:53
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::displacementLaplacianFvMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: displacementLaplacianFvMotionSolver.C:262
motionDiffusivity.H
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
found
bool found
Definition: TABSMDCalcMethod2.H:32
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::points0MotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Definition: points0MotionSolver.C:157
Foam::fvMatrix::solverDict
const dictionary & solverDict() const
Return the solver dictionary taking into account finalIteration.
Definition: fvMatrix.C:1649
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::GeometricField< vector, pointPatchField, pointMesh >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
motionInterpolation.H