mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
INT: Integration of isoAdvector and supporting material
Community contribution from Johan Roenby, DHI
IsoAdvector is a geometric Volume-of-Fluid method for advection of a
sharp interface between two incompressible fluids. It works on both
structured and unstructured meshes with no requirements on cell shapes.
IsoAdvector is as an alternative choice for the interface compression
treatment with the MULES limiter implemented in the interFoam family
of solvers.
The isoAdvector concept and code was developed at DHI and was funded
by a Sapere Aude postdoc grant to Johan Roenby from The Danish Council
for Independent Research | Technology and Production Sciences (Grant-ID:
DFF - 1337-00118B - FTP).
Co-funding is also provided by the GTS grant to DHI from the Danish
Agency for Science, Technology and Innovation.
The ideas behind and performance of the isoAdvector scheme is
documented in:
Roenby J, Bredmose H, Jasak H. 2016 A computational method for sharp
interface advection. R. Soc. open sci. 3: 160405.
[http://dx.doi.org/10.1098/rsos.160405](http://dx.doi.org/10.1098/rsos.160405)
Videos showing isoAdvector's performance with a number of standard
test cases can be found in this youtube channel:
https://www.youtube.com/channel/UCt6Idpv4C8TTgz1iUX0prAA
Project contributors:
* Johan Roenby <jro@dhigroup.com> (Inventor and main developer)
* Hrvoje Jasak <hrvoje.jasak@fsb.hr> (Consistent treatment of
boundary faces including processor boundaries, parallelisation,
code clean up
* Henrik Bredmose <hbre@dtu.dk> (Assisted in the conceptual
development)
* Vuko Vukcevic <vuko.vukcevic@fsb.hr> (Code review, profiling,
porting to foam-extend, bug fixing, testing)
* Tomislav Maric <tomislav@sourceflux.de> (Source file
rearrangement)
* Andy Heather <a.heather@opencfd.co.uk> (Integration into OpenFOAM
for v1706 release)
See the integration repository below to see the full set of changes
implemented for release into OpenFOAM v1706
https://develop.openfoam.com/Community/Integration-isoAdvector
This commit is contained in:
@ -18,6 +18,8 @@ nearWallFields/findCellParticleCloud.C
|
||||
processorField/processorField.C
|
||||
readFields/readFields.C
|
||||
|
||||
setFlow/setFlow.C
|
||||
|
||||
streamLine/streamLine.C
|
||||
streamLine/streamLineBase.C
|
||||
streamLine/streamLineParticle.C
|
||||
|
||||
448
src/functionObjects/field/setFlow/setFlow.C
Normal file
448
src/functionObjects/field/setFlow/setFlow.C
Normal file
@ -0,0 +1,448 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
|
||||
\\/ 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 "setFlow.H"
|
||||
#include "volFields.H"
|
||||
#include "surfaceFields.H"
|
||||
#include "fvcFlux.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "fvcSurfaceIntegrate.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace functionObjects
|
||||
{
|
||||
defineTypeNameAndDebug(setFlow, 0);
|
||||
|
||||
addToRunTimeSelectionTable
|
||||
(
|
||||
functionObject,
|
||||
setFlow,
|
||||
dictionary
|
||||
);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
const Foam::Enum<Foam::functionObjects::setFlow::modeType>
|
||||
Foam::functionObjects::setFlow::modeTypeNames
|
||||
{
|
||||
{ functionObjects::setFlow::modeType::FUNCTION, "function" },
|
||||
{ functionObjects::setFlow::modeType::ROTATION, "rotation" },
|
||||
{ functionObjects::setFlow::modeType::VORTEX2D, "vortex2D" },
|
||||
{ functionObjects::setFlow::modeType::VORTEX3D, "vortex3D" }
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
void Foam::functionObjects::setFlow::setPhi(const volVectorField& U)
|
||||
{
|
||||
surfaceScalarField* phiptr =
|
||||
mesh_.lookupObjectRefPtr<surfaceScalarField>(phiName_);
|
||||
|
||||
if (!phiptr)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
if (rhoName_ != "none")
|
||||
{
|
||||
const volScalarField* rhoptr =
|
||||
mesh_.lookupObjectPtr<volScalarField>(rhoName_);
|
||||
|
||||
if (rhoptr)
|
||||
{
|
||||
const volScalarField& rho = *rhoptr;
|
||||
*phiptr = fvc::flux(rho*U);
|
||||
}
|
||||
else
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Unable to find rho field'" << rhoName_
|
||||
<< "' in the mesh database. Available fields are:"
|
||||
<< mesh_.names<volScalarField>()
|
||||
<< exit(FatalError);
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
*phiptr = fvc::flux(U);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::functionObjects::setFlow::setFlow
|
||||
(
|
||||
const word& name,
|
||||
const Time& runTime,
|
||||
const dictionary& dict
|
||||
)
|
||||
:
|
||||
fvMeshFunctionObject(name, runTime, dict),
|
||||
UName_("U"),
|
||||
rhoName_("none"),
|
||||
phiName_("phi"),
|
||||
mode_(modeType::FUNCTION),
|
||||
reverseTime_(VGREAT),
|
||||
scalePtr_(nullptr),
|
||||
origin_(vector::zero),
|
||||
R_(tensor::I),
|
||||
omegaPtr_(nullptr),
|
||||
velocityPtr_(nullptr)
|
||||
{
|
||||
read(dict);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::functionObjects::setFlow::~setFlow()
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
bool Foam::functionObjects::setFlow::read(const dictionary& dict)
|
||||
{
|
||||
if (fvMeshFunctionObject::read(dict))
|
||||
{
|
||||
Info<< name() << ":" << endl;
|
||||
mode_ = modeTypeNames.read(dict.lookup("mode"));
|
||||
|
||||
Info<< " operating mode: " << modeTypeNames[mode_] << endl;
|
||||
|
||||
if (dict.readIfPresent("U", UName_))
|
||||
{
|
||||
Info<< " U field name: " << UName_ << endl;
|
||||
}
|
||||
|
||||
if (dict.readIfPresent("rho", rhoName_))
|
||||
{
|
||||
Info<< " rho field name: " << rhoName_ << endl;
|
||||
}
|
||||
|
||||
if (dict.readIfPresent("phi", phiName_))
|
||||
{
|
||||
Info<< " phi field name: " << phiName_ << endl;
|
||||
}
|
||||
|
||||
if (dict.readIfPresent("reverseTime", reverseTime_))
|
||||
{
|
||||
Info<< " reverse flow direction at time: " << reverseTime_
|
||||
<< endl;
|
||||
reverseTime_ = mesh_.time().userTimeToTime(reverseTime_);
|
||||
}
|
||||
|
||||
// Scaling is applied across all modes
|
||||
scalePtr_ = Function1<scalar>::New("scale", dict);
|
||||
|
||||
switch (mode_)
|
||||
{
|
||||
case modeType::FUNCTION:
|
||||
{
|
||||
velocityPtr_ = Function1<vector>::New("velocity", dict);
|
||||
break;
|
||||
}
|
||||
case modeType::ROTATION:
|
||||
{
|
||||
omegaPtr_ = Function1<scalar>::New("omega", dict);
|
||||
dict.lookup("origin") >> origin_;
|
||||
vector refDir(dict.lookup("refDir"));
|
||||
refDir /= mag(refDir) + ROOTVSMALL;
|
||||
vector axis(dict.lookup("axis"));
|
||||
axis /= mag(axis) + ROOTVSMALL;
|
||||
R_ = tensor(refDir, axis, refDir^axis);
|
||||
break;
|
||||
}
|
||||
case modeType::VORTEX2D:
|
||||
case modeType::VORTEX3D:
|
||||
{
|
||||
dict.lookup("origin") >> origin_;
|
||||
vector refDir(dict.lookup("refDir"));
|
||||
refDir /= mag(refDir) + ROOTVSMALL;
|
||||
vector axis(dict.lookup("axis"));
|
||||
axis /= mag(axis) + ROOTVSMALL;
|
||||
R_ = tensor(refDir, axis, refDir^axis);
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
Info<< endl;
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
return false;
|
||||
}
|
||||
|
||||
|
||||
bool Foam::functionObjects::setFlow::execute()
|
||||
{
|
||||
volVectorField* Uptr = mesh_.lookupObjectRefPtr<volVectorField>(UName_);
|
||||
|
||||
surfaceScalarField* phiptr =
|
||||
mesh_.lookupObjectRefPtr<surfaceScalarField>(phiName_);
|
||||
|
||||
Log << nl << name() << ":" << nl;
|
||||
|
||||
if (!Uptr || !phiptr)
|
||||
{
|
||||
Info<< " Either field " << UName_ << " or " << phiName_
|
||||
<< " not found in the mesh database" << nl;
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
const scalar t = mesh_.time().timeOutputValue();
|
||||
|
||||
Log << " setting " << UName_ << " and " << phiName_ << nl;
|
||||
|
||||
volVectorField& U = *Uptr;
|
||||
surfaceScalarField& phi = *phiptr;
|
||||
|
||||
switch (mode_)
|
||||
{
|
||||
case modeType::FUNCTION:
|
||||
{
|
||||
const vector Uc = velocityPtr_->value(t);
|
||||
U == dimensionedVector("Uc", dimVelocity, Uc);
|
||||
U.correctBoundaryConditions();
|
||||
setPhi(U);
|
||||
|
||||
break;
|
||||
}
|
||||
case modeType::ROTATION:
|
||||
{
|
||||
const volVectorField& C = mesh_.C();
|
||||
const volVectorField d
|
||||
(
|
||||
typeName + ":d",
|
||||
C - dimensionedVector("origin", dimLength, origin_)
|
||||
);
|
||||
const scalarField x(d.component(vector::X));
|
||||
const scalarField z(d.component(vector::Z));
|
||||
|
||||
const scalar omega = omegaPtr_->value(t);
|
||||
vectorField& Uc = U.primitiveFieldRef();
|
||||
Uc.replace(vector::X, -omega*z);
|
||||
Uc.replace(vector::Y, scalar(0));
|
||||
Uc.replace(vector::Z, omega*x);
|
||||
|
||||
volVectorField::Boundary& Ubf = U.boundaryFieldRef();
|
||||
forAll(Ubf, patchi)
|
||||
{
|
||||
vectorField& Uf = Ubf[patchi];
|
||||
if (Uf.size())
|
||||
{
|
||||
const vectorField& Cp = C.boundaryField()[patchi];
|
||||
const vectorField dp(Cp - origin_);
|
||||
const scalarField xp(dp.component(vector::X));
|
||||
const scalarField zp(dp.component(vector::Z));
|
||||
Uf.replace(vector::X, -omega*zp);
|
||||
Uf.replace(vector::Y, scalar(0));
|
||||
Uf.replace(vector::Z, omega*xp);
|
||||
}
|
||||
}
|
||||
|
||||
U = U & R_;
|
||||
U.correctBoundaryConditions();
|
||||
setPhi(U);
|
||||
|
||||
break;
|
||||
}
|
||||
case modeType::VORTEX2D:
|
||||
{
|
||||
const scalar pi = Foam::constant::mathematical::pi;
|
||||
|
||||
const volVectorField& C = mesh_.C();
|
||||
|
||||
const volVectorField d
|
||||
(
|
||||
typeName + ":d",
|
||||
C - dimensionedVector("origin", dimLength, origin_)
|
||||
);
|
||||
const scalarField x(d.component(vector::X));
|
||||
const scalarField z(d.component(vector::Z));
|
||||
|
||||
vectorField& Uc = U.primitiveFieldRef();
|
||||
Uc.replace(vector::X, -sin(2*pi*z)*sqr(sin(pi*x)));
|
||||
Uc.replace(vector::Y, scalar(0));
|
||||
Uc.replace(vector::Z, sin(2*pi*x)*sqr(sin(pi*z)));
|
||||
|
||||
U = U & R_;
|
||||
|
||||
// Calculating incompressible flux based on stream function
|
||||
// Note: R_ rotation not implemented in phi calculation
|
||||
const scalarField xp(mesh_.points().component(0) - origin_[0]);
|
||||
const scalarField yp(mesh_.points().component(1) - origin_[1]);
|
||||
const scalarField zp(mesh_.points().component(2) - origin_[2]);
|
||||
const scalarField psi((1.0/pi)*sqr(sin(pi*xp))*sqr(sin(pi*zp)));
|
||||
|
||||
scalarField& phic = phi.primitiveFieldRef();
|
||||
forAll(phic, fi)
|
||||
{
|
||||
phic[fi] = 0;
|
||||
const face& f = mesh_.faces()[fi];
|
||||
const label nPoints = f.size();
|
||||
|
||||
forAll(f, fpi)
|
||||
{
|
||||
const label p1 = f[fpi];
|
||||
const label p2 = f[(fpi + 1) % nPoints];
|
||||
phic[fi] += 0.5*(psi[p1] + psi[p2])*(yp[p2] - yp[p1]);
|
||||
}
|
||||
}
|
||||
|
||||
surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
|
||||
forAll(phibf, patchi)
|
||||
{
|
||||
scalarField& phif = phibf[patchi];
|
||||
const label start = mesh_.boundaryMesh()[patchi].start();
|
||||
|
||||
forAll(phif, fi)
|
||||
{
|
||||
phif[fi] = 0;
|
||||
const face& f = mesh_.faces()[start + fi];
|
||||
const label nPoints = f.size();
|
||||
|
||||
forAll(f, fpi)
|
||||
{
|
||||
const label p1 = f[fpi];
|
||||
const label p2 = f[(fpi + 1) % nPoints];
|
||||
phif[fi] += 0.5*(psi[p1] + psi[p2])*(yp[p2] - yp[p1]);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
break;
|
||||
}
|
||||
case modeType::VORTEX3D:
|
||||
{
|
||||
const scalar pi = Foam::constant::mathematical::pi;
|
||||
const volVectorField& C = mesh_.C();
|
||||
|
||||
const volVectorField d
|
||||
(
|
||||
typeName + ":d",
|
||||
C - dimensionedVector("origin", dimLength, origin_)
|
||||
);
|
||||
const scalarField x(d.component(vector::X));
|
||||
const scalarField y(d.component(vector::Y));
|
||||
const scalarField z(d.component(vector::Z));
|
||||
|
||||
vectorField& Uc = U.primitiveFieldRef();
|
||||
Uc.replace(vector::X, 2*sqr(sin(pi*x))*sin(2*pi*y)*sin(2*pi*z));
|
||||
Uc.replace(vector::Y, -sin(2*pi*x)*sqr(sin(pi*y))*sin(2*pi*z));
|
||||
Uc.replace(vector::Z, -sin(2*pi*x)*sin(2*pi*y)*sqr(sin(pi*z)));
|
||||
|
||||
U = U & R_;
|
||||
U.correctBoundaryConditions();
|
||||
|
||||
// Calculating phi
|
||||
// Note: R_ rotation not implemented in phi calculation
|
||||
const vectorField Cf(mesh_.Cf().primitiveField() - origin_);
|
||||
const scalarField Xf(Cf.component(vector::X));
|
||||
const scalarField Yf(Cf.component(vector::Y));
|
||||
const scalarField Zf(Cf.component(vector::Z));
|
||||
vectorField Uf(Xf.size());
|
||||
Uf.replace(0, 2*sqr(sin(pi*Xf))*sin(2*pi*Yf)*sin(2*pi*Zf));
|
||||
Uf.replace(1, -sin(2*pi*Xf)*sqr(sin(pi*Yf))*sin(2*pi*Zf));
|
||||
Uf.replace(2, -sin(2*pi*Xf)*sin(2*pi*Yf)*sqr(sin(pi*Zf)));
|
||||
|
||||
scalarField& phic = phi.primitiveFieldRef();
|
||||
const vectorField& Sfc = mesh_.Sf().primitiveField();
|
||||
phic = Uf & Sfc;
|
||||
|
||||
surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
|
||||
const surfaceVectorField::Boundary& Sfbf =
|
||||
mesh_.Sf().boundaryField();
|
||||
const surfaceVectorField::Boundary& Cfbf =
|
||||
mesh_.Cf().boundaryField();
|
||||
|
||||
forAll(phibf, patchi)
|
||||
{
|
||||
scalarField& phif = phibf[patchi];
|
||||
const vectorField& Sff = Sfbf[patchi];
|
||||
const vectorField& Cff = Cfbf[patchi];
|
||||
const scalarField xf(Cff.component(vector::X));
|
||||
const scalarField yf(Cff.component(vector::Y));
|
||||
const scalarField zf(Cff.component(vector::Z));
|
||||
vectorField Uf(xf.size());
|
||||
Uf.replace(0, 2*sqr(sin(pi*xf))*sin(2*pi*yf)*sin(2*pi*zf));
|
||||
Uf.replace(1, -sin(2*pi*xf)*sqr(sin(pi*yf))*sin(2*pi*zf));
|
||||
Uf.replace(2, -sin(2*pi*xf)*sin(2*pi*yf)*sqr(sin(pi*zf)));
|
||||
phif = Uf & Sff;
|
||||
}
|
||||
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if (t > reverseTime_)
|
||||
{
|
||||
Log << " flow direction: reverse" << nl;
|
||||
U.negate();
|
||||
phi.negate();
|
||||
}
|
||||
|
||||
// Apply scaling
|
||||
const scalar s = scalePtr_->value(t);
|
||||
U *= s;
|
||||
phi *= s;
|
||||
|
||||
U.correctBoundaryConditions();
|
||||
|
||||
const scalarField sumPhi(fvc::surfaceIntegrate(phi));
|
||||
Log << " Continuity error: max(mag(sum(phi))) = "
|
||||
<< gMax(mag(sumPhi)) << nl << endl;
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
bool Foam::functionObjects::setFlow::write()
|
||||
{
|
||||
const auto& Uptr = mesh_.lookupObjectRefPtr<volVectorField>(UName_);
|
||||
if (Uptr)
|
||||
{
|
||||
Uptr->write();
|
||||
}
|
||||
|
||||
const auto& phiptr = mesh_.lookupObjectRefPtr<surfaceScalarField>(phiName_);
|
||||
if (phiptr)
|
||||
{
|
||||
phiptr->write();
|
||||
}
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
199
src/functionObjects/field/setFlow/setFlow.H
Normal file
199
src/functionObjects/field/setFlow/setFlow.H
Normal file
@ -0,0 +1,199 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
|
||||
\\/ 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/>.
|
||||
|
||||
Class
|
||||
Foam::functionObjects::setFlow
|
||||
|
||||
Group
|
||||
grpFieldFunctionObjects
|
||||
|
||||
Description
|
||||
Provides options to set the velocity and flux fields as a function of time
|
||||
|
||||
Useful for testing advection behaviour of numerical schemes by e.g.
|
||||
imposing solid body rotation, vortex flows. All types include a scaling
|
||||
Foam::Function1 type enabling the strength of the transformation to vary as
|
||||
a function of time.
|
||||
|
||||
Usage
|
||||
Example of function object specification:
|
||||
\verbatim
|
||||
setFlow1
|
||||
{
|
||||
type setFlow;
|
||||
libs ("libfieldFunctionObjects.so");
|
||||
...
|
||||
mode rotation;
|
||||
scale 1;
|
||||
reverseTime 1;
|
||||
omega 6.28318530718;
|
||||
origin (0.5 0 0.5);
|
||||
refDir (1 0 0);
|
||||
axis (0 1 0);
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
Where the entries comprise:
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
type | Type name: setFlow | yes |
|
||||
U | Name of velocity field | no | U
|
||||
rho | Name of density field | no | none
|
||||
phi | Name of flux field | no | phi
|
||||
mode | operating mode - see below | yes |
|
||||
\endtable
|
||||
|
||||
Available \c mode types include:
|
||||
- function
|
||||
- rortation
|
||||
- vortex2D
|
||||
- vortex3D
|
||||
|
||||
|
||||
See also
|
||||
Foam::functionObjects::fvMeshFunctionObject
|
||||
|
||||
SourceFiles
|
||||
setFlow.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef functionObjects_setFlow_H
|
||||
#define functionObjects_setFlow_H
|
||||
|
||||
#include "fvMeshFunctionObject.H"
|
||||
#include "Function1.H"
|
||||
#include "Enum.H"
|
||||
#include "point.H"
|
||||
#include "volFieldsFwd.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace functionObjects
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class setFlow Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class setFlow
|
||||
:
|
||||
public fvMeshFunctionObject
|
||||
{
|
||||
enum class modeType
|
||||
{
|
||||
FUNCTION,
|
||||
ROTATION,
|
||||
VORTEX2D,
|
||||
VORTEX3D
|
||||
};
|
||||
|
||||
static const Foam::Enum<modeType> modeTypeNames;
|
||||
|
||||
|
||||
// Private Data
|
||||
|
||||
//- Name of vecloity field, default = "U"
|
||||
word UName_;
|
||||
|
||||
//- Name of density field, default = "none"
|
||||
word rhoName_;
|
||||
|
||||
//- Name of flux field, default = "phi"
|
||||
word phiName_;
|
||||
|
||||
//- Operating mode
|
||||
modeType mode_;
|
||||
|
||||
//- Reverse time
|
||||
scalar reverseTime_;
|
||||
|
||||
//- Scaling function
|
||||
autoPtr<Function1<scalar>> scalePtr_;
|
||||
|
||||
|
||||
// Rotation
|
||||
|
||||
//- Origin
|
||||
point origin_;
|
||||
|
||||
//- Rotation tensor for rotational mode
|
||||
tensor R_;
|
||||
|
||||
//- Rotational speed function
|
||||
autoPtr<Function1<scalar>> omegaPtr_;
|
||||
|
||||
|
||||
// Function
|
||||
|
||||
//- Velocity function
|
||||
autoPtr<Function1<vector>> velocityPtr_;
|
||||
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
//- Set the flux field
|
||||
void setPhi(const volVectorField& U);
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("setFlow");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from Time and dictionary
|
||||
setFlow
|
||||
(
|
||||
const word& name,
|
||||
const Time& runTime,
|
||||
const dictionary& dict
|
||||
);
|
||||
|
||||
|
||||
//- Destructor
|
||||
virtual ~setFlow();
|
||||
|
||||
|
||||
virtual bool read(const dictionary& dict);
|
||||
|
||||
virtual bool execute();
|
||||
|
||||
virtual bool write();
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace functionObjects
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
Reference in New Issue
Block a user