velocityComponentLaplacianFvMotionSolver.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 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(velocityComponentLaplacianFvMotionSolver, 0);
41 
43  (
44  motionSolver,
45  velocityComponentLaplacianFvMotionSolver,
46  dictionary
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
53 Foam::velocityComponentLaplacianFvMotionSolver::
54 velocityComponentLaplacianFvMotionSolver
55 (
56  const polyMesh& mesh,
57  const IOdictionary& dict
58 )
59 :
62  cellMotionU_
63  (
64  IOobject
65  (
66  "cellMotionU" + cmptName_,
67  mesh.time().timeName(),
68  mesh,
69  IOobject::READ_IF_PRESENT,
70  IOobject::AUTO_WRITE
71  ),
72  fvMesh_,
73  dimensionedScalar(pointMotionU_.dimensions(), Zero),
74  cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
75  ),
76  interpolationPtr_
77  (
78  coeffDict().found("interpolation")
79  ? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
80  : motionInterpolation::New(fvMesh_)
81  ),
82  diffusivityPtr_
83  (
84  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
85  )
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
100 {
101  interpolationPtr_->interpolate
102  (
103  cellMotionU_,
104  pointMotionU_
105  );
106 
107  tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
108 
109  tcurPoints.ref().replace
110  (
111  cmpt_,
112  tcurPoints().component(cmpt_)
113  + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
114  );
116  twoDCorrectPoints(tcurPoints.ref());
117 
118  return tcurPoints;
119 }
120 
121 
123 {
124  // The points have moved so before interpolation update
125  // the fvMotionSolver accordingly
126  movePoints(fvMesh_.points());
127 
128  diffusivityPtr_->correct();
129  pointMotionU_.boundaryFieldRef().updateCoeffs();
130 
132 
133  // We explicitly do NOT want to interpolate the motion inbetween
134  // different regions so bypass all the matrix manipulation.
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 
155 (
156  const mapPolyMesh& mpm
157 )
158 {
160 
161  // Update diffusivity. Note two stage to make sure old one is de-registered
162  // before creating/registering new one.
163  diffusivityPtr_.reset(nullptr);
164  diffusivityPtr_ = motionDiffusivity::New
165  (
166  fvMesh_,
167  coeffDict().lookup("diffusivity")
168  );
169 }
170 
171 
172 // ************************************************************************* //
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
dictionary dict
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 tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Finite-volume options.
Definition: fvOptions.H:51
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
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
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Virtual base class for velocity motion solver.
dynamicFvMesh & mesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
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.
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)
virtual void updateMesh(const mapPolyMesh &)
Update topology.
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
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))
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
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
bool found
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127