velocityLaplacianFvMotionSolver.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-2016 OpenFOAM Foundation
9  Copyright (C) 2016-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 "fvOptions.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(velocityLaplacianFvMotionSolver, 0);
41 
43  (
44  motionSolver,
45  velocityLaplacianFvMotionSolver,
46  dictionary
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
53 Foam::velocityLaplacianFvMotionSolver::velocityLaplacianFvMotionSolver
54 (
55  const polyMesh& mesh,
56  const IOdictionary& dict
57 )
58 :
59  velocityMotionSolver(mesh, dict, typeName),
61  cellMotionU_
62  (
63  IOobject
64  (
65  "cellMotionU",
66  mesh.time().timeName(),
67  mesh,
68  IOobject::READ_IF_PRESENT,
69  IOobject::AUTO_WRITE
70  ),
71  fvMesh_,
72  dimensionedVector(pointMotionU_.dimensions(), Zero),
73  cellMotionBoundaryTypes<vector>(pointMotionU_.boundaryField())
74  ),
75  interpolationPtr_
76  (
77  coeffDict().found("interpolation")
78  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
79  : motionInterpolation::New(fvMesh_)
80  ),
81  diffusivityPtr_
82  (
83  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
84  )
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
98 {
99  interpolationPtr_->interpolate
100  (
101  cellMotionU_,
102  pointMotionU_
103  );
104 
105  tmp<pointField> tcurPoints
106  (
107  fvMesh_.points()
108  + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
109  );
111  twoDCorrectPoints(tcurPoints.ref());
112 
113  return tcurPoints;
114 }
115 
116 
118 {
119  // The points have moved so before interpolation update
120  // the fvMotionSolver accordingly
121  movePoints(fvMesh_.points());
122 
123  diffusivityPtr_->correct();
124  pointMotionU_.boundaryFieldRef().updateCoeffs();
125 
127 
128  const label nNonOrthCorr
129  (
130  getOrDefault<label>("nNonOrthogonalCorrectors", 1)
131  );
132 
133  for (label i=0; i<nNonOrthCorr; ++i)
134  {
136  (
138  (
139  dimensionedScalar("viscosity", dimViscosity, 1.0)
140  * diffusivityPtr_->operator()(),
141  cellMotionU_,
142  "laplacian(diffusivity,cellMotionU)"
143  )
144  ==
145  fvOptions(cellMotionU_)
146  );
147 
150  fvOptions.correct(cellMotionU_);
151  }
152 }
153 
154 
156 (
157  const mapPolyMesh& mpm
158 )
159 {
161 
162  // Update diffusivity. Note two stage to make sure old one is de-registered
163  // before creating/registering new one.
164  diffusivityPtr_.reset(nullptr);
165  diffusivityPtr_ = motionDiffusivity::New
166  (
167  fvMesh_,
168  coeffDict().lookup("diffusivity")
169  );
170 }
171 
172 
173 // ************************************************************************* //
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
dictionary dict
Virtual base class for velocity motion solver.
const dimensionSet dimViscosity
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 void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Finite-volume options.
Definition: fvOptions.H:51
Lookup type of boundary radiation properties.
Definition: lookup.H:57
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
const int nNonOrthCorr
SolverPerformance< Type > solveSegregatedOrCoupled(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:55
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
Base class for fvMesh based motionSolvers.
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
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
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...
Abstract base class for cell-centre mesh motion diffusivity.
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)
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
fvVectorMatrix & UEqn
Definition: UEqn.H:13
virtual void updateMesh(const mapPolyMesh &)
Update topology.
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.
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