Files
openfoam/src/fvMotionSolver/fvMotionSolvers/componentVelocity/componentLaplacian/velocityComponentLaplacianFvMotionSolver.C
2012-08-31 16:54:55 +01:00

157 lines
4.2 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "velocityComponentLaplacianFvMotionSolver.H"
#include "motionDiffusivity.H"
#include "fvmLaplacian.H"
#include "addToRunTimeSelectionTable.H"
#include "volPointInterpolation.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(velocityComponentLaplacianFvMotionSolver, 0);
addToRunTimeSelectionTable
(
motionSolver,
velocityComponentLaplacianFvMotionSolver,
dictionary
);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::velocityComponentLaplacianFvMotionSolver::
velocityComponentLaplacianFvMotionSolver
(
const polyMesh& mesh,
const IOdictionary& dict
)
:
componentVelocityMotionSolver(mesh, dict, typeName),
fvMotionSolverCore(mesh),
cellMotionU_
(
IOobject
(
"cellMotionU" + cmptName_,
mesh.time().timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
fvMesh_,
dimensionedScalar
(
"cellMotionU",
pointMotionU_.dimensions(),
0
),
cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
),
diffusivityPtr_
(
motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::velocityComponentLaplacianFvMotionSolver::
~velocityComponentLaplacianFvMotionSolver()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::pointField>
Foam::velocityComponentLaplacianFvMotionSolver::curPoints() const
{
volPointInterpolation::New(fvMesh_).interpolate
(
cellMotionU_,
pointMotionU_
);
tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
tcurPoints().replace
(
cmpt_,
tcurPoints().component(cmpt_)
+ fvMesh_.time().deltaTValue()*pointMotionU_.internalField()
);
twoDCorrectPoints(tcurPoints());
return tcurPoints;
}
void Foam::velocityComponentLaplacianFvMotionSolver::solve()
{
// The points have moved so before interpolation update
// the fvMotionSolver accordingly
movePoints(fvMesh_.points());
diffusivityPtr_->correct();
pointMotionU_.boundaryField().updateCoeffs();
Foam::solve
(
fvm::laplacian
(
diffusivityPtr_->operator()(),
cellMotionU_,
"laplacian(diffusivity,cellMotionU)"
)
);
}
void Foam::velocityComponentLaplacianFvMotionSolver::updateMesh
(
const mapPolyMesh& mpm
)
{
componentVelocityMotionSolver::updateMesh(mpm);
// Update diffusivity. Note two stage to make sure old one is de-registered
// before creating/registering new one.
diffusivityPtr_.reset(NULL);
diffusivityPtr_ = motionDiffusivity::New
(
fvMesh_,
coeffDict().lookup("diffusivity")
);
}
// ************************************************************************* //