Updates to MRF handling for block-coupled solver

This commit is contained in:
Hrvoje Jasak 2017-11-18 12:21:56 +00:00
parent 111979afe5
commit 0413c44e09
7 changed files with 117 additions and 15 deletions

View file

@ -9,13 +9,10 @@
+ turbulence->divDevReff() + turbulence->divDevReff()
); );
// Add MRF sources // Add MRF and porous sources implicitly. HJ, 18/Nov/2017
mrfZones.addCoriolis(UEqn);
// Add porous sources
tmp<volTensorField> tTU; tmp<volTensorField> tTU;
if (addPorosity) if (addMRF || addPorosity)
{ {
tTU = tmp<volTensorField> tTU = tmp<volTensorField>
( (
@ -35,6 +32,10 @@
); );
volTensorField& TU = tTU(); volTensorField& TU = tTU();
// Add implicit MRF source as a Hodge dual of the rotational velocity
TU += *mrfZones.omega();
// Add implicit resistance
pZones.addResistance(UEqn, TU); pZones.addResistance(UEqn, TU);
trTU = inv(TU + tensor(I)*UEqn.A()); trTU = inv(TU + tensor(I)*UEqn.A());
@ -52,7 +53,7 @@
// Insert momentum equation // Insert momentum equation
UpEqn.insertEquation(0, UEqn); UpEqn.insertEquation(0, UEqn);
if (addPorosity) if (addMRF || addPorosity)
{ {
// Manually over-ride the 3x3 block to handle the off-diagonal // Manually over-ride the 3x3 block to handle the off-diagonal
// part of the Ap coefficient // part of the Ap coefficient
@ -62,11 +63,43 @@
const scalarField& V = mesh.V().field(); const scalarField& V = mesh.V().field();
// Note: insertion should only happen in porous cell zones // Note: insertion should only happen in MRF and porous cell zones
// HJ, 14/Mar/2016 // HJ, 14/Mar/2016
// Warning. Possible problem with a zone which is both MRF and porous
// The solution is to do the loop below everywhere
// Currently only inserting zones for efficiency. HJ, 18/Nov/2017
register label cellI; register label cellI;
forAll (mrfZones, mrfZoneI)
{
const labelList& curZoneCells = mrfZones[mrfZoneI].zone();
// Loop over all cells in the zone
forAll (curZoneCells, zcI)
{
cellI = curZoneCells[zcI];
const scalar& cellV = V[cellI];
const tensor& cellTU = TUIn[cellI];
CoeffField<vector4>::squareType& cellDD = DD[cellI];
cellDD(0, 0) += cellV*cellTU.xx();
cellDD(0, 1) += cellV*cellTU.xy();
cellDD(0, 2) += cellV*cellTU.xz();
cellDD(1, 0) += cellV*cellTU.yx();
cellDD(1, 1) += cellV*cellTU.yy();
cellDD(2, 2) += cellV*cellTU.yz();
cellDD(2, 0) += cellV*cellTU.zx();
cellDD(2, 1) += cellV*cellTU.zy();
cellDD(2, 2) += cellV*cellTU.zz();
}
}
forAll (pZones, pZoneI) forAll (pZones, pZoneI)
{ {
const labelList& curZoneCells = pZones[pZoneI].zone(); const labelList& curZoneCells = pZones[pZoneI].zone();

View file

@ -1,2 +1,9 @@
MRFZones mrfZones(mesh); MRFZones mrfZones(mesh);
mrfZones.correctBoundaryVelocity(U); mrfZones.correctBoundaryVelocity(U);
bool addMRF = false;
if (!mrfZones.empty())
{
addMRF = true;
}

View file

@ -2,7 +2,7 @@
tmp<fvScalarMatrix> tpEqn; tmp<fvScalarMatrix> tpEqn;
tmp<surfaceScalarField> tpresSource; tmp<surfaceScalarField> tpresSource;
if (addPorosity) if (addMRF || addPorosity)
{ {
// Collect pressure source with tensorial 1/Ap // Collect pressure source with tensorial 1/Ap
const volTensorField& rTU = trTU(); const volTensorField& rTU = trTU();

View file

@ -581,6 +581,20 @@ void Foam::MRFZone::addCoriolis
} }
void Foam::MRFZone::addOmega(volVectorField& omg) const
{
const vector rotVel = Omega();
// Set omega in all cells of the rotating cell zone
const labelList& cells = mesh_.cellZones()[cellZoneID_];
forAll (cells, i)
{
omg[cells[i]] = rotVel;
}
}
void Foam::MRFZone::relativeVelocity(volVectorField& U) const void Foam::MRFZone::relativeVelocity(volVectorField& U) const
{ {
const volVectorField& C = mesh_.C(); const volVectorField& C = mesh_.C();

View file

@ -203,6 +203,18 @@ public:
return name_; return name_;
} }
//- Return cellZone number
label zoneId() const
{
return cellZoneID_;
}
//- Return cellZone
const cellZone& zone() const
{
return mesh_.cellZones()[cellZoneID_];
}
//- Return the MRFZone dictionary //- Return the MRFZone dictionary
const dictionary& dict() const const dictionary& dict() const
{ {
@ -250,6 +262,9 @@ public:
fvVectorMatrix& UEqn fvVectorMatrix& UEqn
) const; ) const;
//- Add omega within the MRF region
void addOmega(volVectorField& omg) const;
//- Make the given absolute velocity relative within the MRF region //- Make the given absolute velocity relative within the MRF region
void relativeVelocity(volVectorField& U) const; void relativeVelocity(volVectorField& U) const;
@ -276,7 +291,7 @@ public:
surfaceScalarField& phi surfaceScalarField& phi
) const; ) const;
//- Compute the pseudo mesh phi of the MRF region //- Compute pseudo mesh phi for the MRF region
void meshPhi(surfaceScalarField& phi) const; void meshPhi(surfaceScalarField& phi) const;
//- Correct the boundary velocity for the roation of the MRF region //- Correct the boundary velocity for the roation of the MRF region

View file

@ -57,6 +57,35 @@ Foam::MRFZones::MRFZones(const fvMesh& mesh)
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::volVectorField> Foam::MRFZones::omega() const
{
tmp<volVectorField> tMRFZonesOmega
(
new volVectorField
(
IOobject
(
"MRFZonesOmega",
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh_,
dimensionedVector("zero", dimless/dimTime, vector::zero)
)
);
volVectorField& MRFZonesOmega = tMRFZonesOmega();
forAll (*this, i)
{
operator[](i).addOmega(MRFZonesOmega);
}
return tMRFZonesOmega;
}
Foam::tmp<Foam::surfaceScalarField> Foam::MRFZones::fluxCorrection() const Foam::tmp<Foam::surfaceScalarField> Foam::MRFZones::fluxCorrection() const
{ {
tmp<surfaceScalarField> tMRFZonesPhiCorr tmp<surfaceScalarField> tMRFZonesPhiCorr
@ -88,7 +117,7 @@ Foam::tmp<Foam::surfaceScalarField> Foam::MRFZones::fluxCorrection() const
Foam::tmp<Foam::surfaceScalarField> Foam::MRFZones::meshPhi() const Foam::tmp<Foam::surfaceScalarField> Foam::MRFZones::meshPhi() const
{ {
tmp<surfaceScalarField> tMRFZonesPhiCorr tmp<surfaceScalarField> tMRFZonesFaceU
( (
new surfaceScalarField new surfaceScalarField
( (
@ -104,14 +133,14 @@ Foam::tmp<Foam::surfaceScalarField> Foam::MRFZones::meshPhi() const
dimensionedScalar("zero", dimVolume/dimTime, 0) dimensionedScalar("zero", dimVolume/dimTime, 0)
) )
); );
surfaceScalarField& MRFZonesPhiCorr = tMRFZonesPhiCorr(); surfaceScalarField& MRFZonesFaceU = tMRFZonesFaceU();
forAll (*this, i) forAll (*this, i)
{ {
operator[](i).meshPhi(MRFZonesPhiCorr); operator[](i).meshPhi(MRFZonesFaceU);
} }
return tMRFZonesPhiCorr; return tMRFZonesFaceU;
} }

View file

@ -77,12 +77,16 @@ public:
// Member Functions // Member Functions
//- Return cell centre omega over all zones
tmp<volVectorField> omega() const;
//- Return raw correction flux //- Return raw correction flux
tmp<surfaceScalarField> fluxCorrection() const; tmp<surfaceScalarField> fluxCorrection() const;
//- Return pseudo mesh flux //- Return pseudo mesh flux
tmp<surfaceScalarField> meshPhi() const; tmp<surfaceScalarField> meshPhi() const;
// Incompressible MRF // Incompressible MRF
//- Add the Coriolis force contribution to the momentum equation //- Add the Coriolis force contribution to the momentum equation