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 | Copyright (C) 2011-2016 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 
49 Foam::velocityComponentLaplacianFvMotionSolver::
50 velocityComponentLaplacianFvMotionSolver
51 (
52  const polyMesh& mesh,
53  const IOdictionary& dict
54 )
55 :
56  componentVelocityMotionSolver(mesh, dict, typeName),
57  fvMotionSolverCore(mesh),
58  cellMotionU_
59  (
60  IOobject
61  (
62  "cellMotionU" + cmptName_,
63  mesh.time().timeName(),
64  mesh,
67  ),
68  fvMesh_,
70  (
71  "cellMotionU",
72  pointMotionU_.dimensions(),
73  0
74  ),
75  cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
76  ),
77  diffusivityPtr_
78  (
79  motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
80  )
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
85 
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
95 {
97  (
98  cellMotionU_,
99  pointMotionU_
100  );
101 
102  tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
103 
104  tcurPoints.ref().replace
105  (
106  cmpt_,
107  tcurPoints().component(cmpt_)
108  + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
109  );
110 
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  (
129  (
130  diffusivityPtr_->operator()(),
131  cellMotionU_,
132  "laplacian(diffusivity,cellMotionU)"
133  )
134  );
135 }
136 
137 
139 (
140  const mapPolyMesh& mpm
141 )
142 {
144 
145  // Update diffusivity. Note two stage to make sure old one is de-registered
146  // before creating/registering new one.
147  diffusivityPtr_.reset(NULL);
148  diffusivityPtr_ = motionDiffusivity::New
149  (
150  fvMesh_,
151  coeffDict().lookup("diffusivity")
152  );
153 }
154 
155 
156 // ************************************************************************* //
const Time & time() const
Return time.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Calculate the matrix for the laplacian of the field.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:662
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)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
Virtual base class for velocity motion solver.
dynamicFvMesh & mesh
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Base class for fvMesh based motionSolvers.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual void updateMesh(const mapPolyMesh &)
Update topology.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:54
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.