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-2019 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 {
36  defineTypeNameAndDebug(velocityComponentLaplacianFvMotionSolver, 0);
37 
39  (
40  motionSolver,
41  velocityComponentLaplacianFvMotionSolver,
42  dictionary
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
51 (
52  const polyMesh& mesh,
53  const dictionary& dict
54 )
55 :
56  componentVelocityMotionSolver(mesh, dict, typeName),
57  fvMotionSolver(mesh),
58  cellMotionU_
59  (
60  IOobject
61  (
62  "cellMotionU" + cmptName_,
63  mesh.time().timeName(),
64  mesh,
67  ),
68  fvMesh_,
69  dimensionedScalar(pointMotionU_.dimensions(), 0),
70  cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
71  ),
72  diffusivityPtr_
73  (
74  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
75  )
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
80 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
90 {
92  (
93  cellMotionU_,
94  pointMotionU_
95  );
96 
97  tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
98 
99  tcurPoints.ref().replace
100  (
101  cmpt_,
102  tcurPoints().component(cmpt_)
103  + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
104  );
105 
106  twoDCorrectPoints(tcurPoints.ref());
107 
108  return tcurPoints;
109 }
110 
111 
113 {
114  // The points have moved so before interpolation update
115  // the fvMotionSolver accordingly
116  movePoints(fvMesh_.points());
117 
118  diffusivityPtr_->correct();
119  pointMotionU_.boundaryFieldRef().updateCoeffs();
120 
122  (
124  (
125  diffusivityPtr_->operator()(),
126  cellMotionU_,
127  "laplacian(diffusivity,cellMotionU)"
128  )
129  );
130 }
131 
132 
134 (
135  const mapPolyMesh& mpm
136 )
137 {
139 
140  // Update diffusivity. Note two stage to make sure old one is de-registered
141  // before creating/registering new one.
142  diffusivityPtr_.reset(nullptr);
143  diffusivityPtr_ = motionDiffusivity::New
144  (
145  fvMesh_,
146  coeffDict().lookup("diffusivity")
147  );
148 }
149 
150 
151 // ************************************************************************* //
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
velocityComponentLaplacianFvMotionSolver(const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Calculate the matrix for the laplacian of the field.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:472
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
static const volPointInterpolation & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Base class for fvMesh based motionSolvers.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
Virtual base class for velocity motion solver.
dynamicFvMesh & mesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Time & time() const
Return time.
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:46
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.