OLD: pEqn.solve(mesh.solver(p.select(piso.finalInnerIter()))); pEqn.solve(mesh.solver("Yi")); NEW: pEqn.solve(p.select(piso.finalInnerIter())); pEqn.solve("Yi");
173 lines
4.8 KiB
C
173 lines
4.8 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | www.openfoam.com
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
Copyright (C) 2011-2016 OpenFOAM Foundation
|
|
Copyright (C) 2016 OpenCFD Ltd.
|
|
-------------------------------------------------------------------------------
|
|
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 "motionInterpolation.H"
|
|
#include "motionDiffusivity.H"
|
|
#include "fvmLaplacian.H"
|
|
#include "addToRunTimeSelectionTable.H"
|
|
#include "fvOptions.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),
|
|
fvMotionSolver(mesh),
|
|
cellMotionU_
|
|
(
|
|
IOobject
|
|
(
|
|
"cellMotionU" + cmptName_,
|
|
mesh.time().timeName(),
|
|
mesh,
|
|
IOobject::READ_IF_PRESENT,
|
|
IOobject::AUTO_WRITE
|
|
),
|
|
fvMesh_,
|
|
dimensionedScalar(pointMotionU_.dimensions(), Zero),
|
|
cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
|
|
),
|
|
interpolationPtr_
|
|
(
|
|
coeffDict().found("interpolation")
|
|
? motionInterpolation::New(fvMesh_, coeffDict().lookup("interpolation"))
|
|
: motionInterpolation::New(fvMesh_)
|
|
),
|
|
diffusivityPtr_
|
|
(
|
|
motionDiffusivity::New(fvMesh_, coeffDict().lookup("diffusivity"))
|
|
)
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
|
|
Foam::velocityComponentLaplacianFvMotionSolver::
|
|
~velocityComponentLaplacianFvMotionSolver()
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
Foam::tmp<Foam::pointField>
|
|
Foam::velocityComponentLaplacianFvMotionSolver::curPoints() const
|
|
{
|
|
interpolationPtr_->interpolate
|
|
(
|
|
cellMotionU_,
|
|
pointMotionU_
|
|
);
|
|
|
|
tmp<pointField> tcurPoints(new pointField(fvMesh_.points()));
|
|
|
|
tcurPoints.ref().replace
|
|
(
|
|
cmpt_,
|
|
tcurPoints().component(cmpt_)
|
|
+ fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
|
|
);
|
|
|
|
twoDCorrectPoints(tcurPoints.ref());
|
|
|
|
return tcurPoints;
|
|
}
|
|
|
|
|
|
void Foam::velocityComponentLaplacianFvMotionSolver::solve()
|
|
{
|
|
// The points have moved so before interpolation update
|
|
// the fvMotionSolver accordingly
|
|
movePoints(fvMesh_.points());
|
|
|
|
diffusivityPtr_->correct();
|
|
pointMotionU_.boundaryFieldRef().updateCoeffs();
|
|
|
|
fv::options& fvOptions(fv::options::New(fvMesh_));
|
|
|
|
// We explicitly do NOT want to interpolate the motion inbetween
|
|
// different regions so bypass all the matrix manipulation.
|
|
fvScalarMatrix TEqn
|
|
(
|
|
fvm::laplacian
|
|
(
|
|
dimensionedScalar("viscosity", dimViscosity, 1.0)
|
|
*diffusivityPtr_->operator()(),
|
|
cellMotionU_,
|
|
"laplacian(diffusivity,cellMotionU)"
|
|
)
|
|
==
|
|
fvOptions(cellMotionU_)
|
|
);
|
|
|
|
fvOptions.constrain(TEqn);
|
|
TEqn.solveSegregatedOrCoupled();
|
|
fvOptions.correct(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(nullptr);
|
|
diffusivityPtr_ = motionDiffusivity::New
|
|
(
|
|
fvMesh_,
|
|
coeffDict().lookup("diffusivity")
|
|
);
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|