velocityComponentLaplacianFvMotionSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "motionDiffusivity.H"
28 #include "fvmLaplacian.H"
30 #include "volPointInterpolation.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 
39  (
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
51 (
52  const word& name,
53  const polyMesh& mesh,
54  const dictionary& dict
55 )
56 :
59  cellMotionU_
60  (
61  IOobject
62  (
63  "cellMotionU" + cmptName_,
64  mesh.time().name(),
65  mesh,
66  IOobject::READ_IF_PRESENT,
67  IOobject::AUTO_WRITE
68  ),
69  fvMesh_,
70  dimensionedScalar(pointMotionU_.dimensions(), 0),
71  cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
72  ),
73  diffusivityType_(dict.lookup("diffusivity")),
74  diffusivityPtr_
75  (
76  motionDiffusivity::New(fvMesh_, diffusivityType_)
77  )
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
82 
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
92 {
94  (
95  cellMotionU_,
96  pointMotionU_
97  );
98 
99  tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
100 
101  tcurPoints.ref().replace
102  (
103  cmpt_,
104  tcurPoints().component(cmpt_)
105  + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
106  );
107 
108  twoDCorrectPoints(tcurPoints.ref());
109 
110  return tcurPoints;
111 }
112 
113 
115 {
116  // The points have moved so before interpolation update
117  // the fvMotionSolver accordingly
118  movePoints(fvMesh_.points());
119 
120  diffusivityPtr_->correct();
121  pointMotionU_.boundaryFieldRef().updateCoeffs();
122 
124  (
126  (
127  diffusivityPtr_->operator()(),
128  cellMotionU_,
129  "laplacian(diffusivity,cellMotionU)"
130  )
131  );
132 }
133 
134 
136 (
137  const polyTopoChangeMap& map
138 )
139 {
141 
142  // Update diffusivity. Note two stage to make sure old one is de-registered
143  // before creating/registering new one.
144  diffusivityPtr_.reset(nullptr);
145  diffusivityType_.rewind();
146  diffusivityPtr_ = motionDiffusivity::New
147  (
148  fvMesh_,
149  diffusivityType_
150  );
151 }
152 
153 
155 (
156  const polyMeshMap& map
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  diffusivityType_.rewind();
165  diffusivityPtr_ = motionDiffusivity::New
166  (
167  fvMesh_,
168  diffusivityType_
169  );
170 }
171 
172 
173 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static volPointInterpolation & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Virtual base class for velocity motion solver.
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base class for fvMesh based motionSolvers.
Abstract base class for cell-centre mesh motion diffusivity.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Mesh motion solver for an fvMesh. Based on solving the cell-centre Laplacian for the given component ...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
velocityComponentLaplacianFvMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the matrix for the laplacian of the field.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dictionary dict