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-2022 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 :
57  componentVelocityMotionSolver(name, mesh, dict, typeName),
58  fvMotionSolver(mesh),
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  diffusivityPtr_
74  (
75  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
76  )
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
91 {
93  (
94  cellMotionU_,
95  pointMotionU_
96  );
97 
98  tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
99 
100  tcurPoints.ref().replace
101  (
102  cmpt_,
103  tcurPoints().component(cmpt_)
104  + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
105  );
106 
107  twoDCorrectPoints(tcurPoints.ref());
108 
109  return tcurPoints;
110 }
111 
112 
114 {
115  // The points have moved so before interpolation update
116  // the fvMotionSolver accordingly
117  movePoints(fvMesh_.points());
118 
119  diffusivityPtr_->correct();
120  pointMotionU_.boundaryFieldRef().updateCoeffs();
121 
123  (
125  (
126  diffusivityPtr_->operator()(),
127  cellMotionU_,
128  "laplacian(diffusivity,cellMotionU)"
129  )
130  );
131 }
132 
133 
135 (
136  const polyTopoChangeMap& map
137 )
138 {
140 
141  // Update diffusivity. Note two stage to make sure old one is de-registered
142  // before creating/registering new one.
143  diffusivityPtr_.reset(nullptr);
144  diffusivityPtr_ = motionDiffusivity::New
145  (
146  fvMesh_,
147  coeffDict().lookup("diffusivity")
148  );
149 }
150 
151 
153 (
154  const polyMeshMap& map
155 )
156 {
158 
159  // Update diffusivity. Note two stage to make sure old one is de-registered
160  // before creating/registering new one.
161  diffusivityPtr_.reset(nullptr);
162  diffusivityPtr_ = motionDiffusivity::New
163  (
164  fvMesh_,
165  coeffDict().lookup("diffusivity")
166  );
167 }
168 
169 
170 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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:181
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
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.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dictionary dict