Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

calculate flux of rhoY in dfHighSpeedFoam #517

Merged
merged 3 commits into from
Oct 12, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
39 changes: 37 additions & 2 deletions applications/solvers/dfHighSpeedFoam/createFields.H
Original file line number Diff line number Diff line change
Expand Up @@ -124,7 +124,20 @@ volVectorField rhoU
rho*U
);

surfaceScalarField phi("phi", fvc::flux(rhoU));
// surfaceScalarField phi("phi", fvc::flux(rhoU));
surfaceScalarField phi
(
IOobject
(
"phi",
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionedScalar("0", dimVelocity*dimArea, 0.0)
);

Info<< "Creating flux properties\n" << endl;
surfaceScalarField rhoPhi
Expand Down Expand Up @@ -271,6 +284,28 @@ forAll(rhoYi,i)
);
}

PtrList<surfaceScalarField> rhoPhiYi(nspecies);
forAll(rhoPhiYi,i)
{
rhoPhiYi.set
(
i,
new surfaceScalarField
(
IOobject
(
"rhoPhi" + Y[i].name(),
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionedScalar("0", dimDensity*dimVelocity*dimArea, 0.0)
)
);
}

#include "createMRF.H"
#include "createClouds.H"

Expand Down Expand Up @@ -339,4 +374,4 @@ volScalarField pMax
IOobject::AUTO_WRITE
),
thermo.p()
);
);
4 changes: 2 additions & 2 deletions applications/solvers/dfHighSpeedFoam/dfHighSpeedFoam.C
Original file line number Diff line number Diff line change
Expand Up @@ -170,7 +170,7 @@ int main(int argc, char *argv[])
volScalarField rPsi("rPsi", 1.0/psi);
volScalarField c("c", sqrt(thermo.Cp()/thermo.Cv()*rPsi));

fluxSchemeFields->update(rho,U,ea,p,c,phi,rhoPhi,rhoUPhi,rhoEPhi);
fluxSchemeFields->update(rho,rhoYi,nspecies,U,ea,p,c,phi,rhoPhi,rhoPhiYi,rhoUPhi,rhoEPhi);

volScalarField muEff("muEff", turbulence->muEff());
volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
Expand Down Expand Up @@ -222,7 +222,7 @@ int main(int argc, char *argv[])
volScalarField rPsi("rPsi", 1.0/psi);
volScalarField c("c", sqrt(thermo.Cp()/thermo.Cv()*rPsi));

fluxSchemeFields->update(rho,U,ea,p,c,phi,rhoPhi,rhoUPhi,rhoEPhi);
fluxSchemeFields->update(rho,rhoYi,nspecies,U,ea,p,c,phi,rhoPhi,rhoPhiYi,rhoUPhi,rhoEPhi);

volScalarField muEff("muEff", turbulence->muEff());
volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
Expand Down
82 changes: 41 additions & 41 deletions applications/solvers/dfHighSpeedFoam/rhoYEqn.H
Original file line number Diff line number Diff line change
Expand Up @@ -53,7 +53,7 @@ time_monitor_corrDiff += double(end - start) / double(CLOCKS_PER_SEC);
forAll(Y, i)
{
volScalarField& Yi = Y[i];
surfaceScalarField rhoPhiYi(rhoPhi*fvc::interpolate(Yi, "Yi").ref());
// surfaceScalarField rhoPhiYi(rhoPhi*fvc::interpolate(Yi, "Yi").ref());

if (!inviscid)
{
Expand All @@ -65,7 +65,7 @@ time_monitor_corrDiff += double(end - start) / double(CLOCKS_PER_SEC);
{
if((ddtSchemes == "RK2SSP") || (ddtSchemes == "RK3SSP"))
{
rhoYi_rhs = -fvc::div(rhoPhiYi);
rhoYi_rhs = -fvc::div(rhoPhiYi[i]);

if(chemScheme == "direct")
{
Expand Down Expand Up @@ -107,62 +107,62 @@ time_monitor_corrDiff += double(end - start) / double(CLOCKS_PER_SEC);
else
{
// original convection term
fvScalarMatrix YiEqn
(
fvm::ddt(rho, Yi)
+ mvConvection->fvmDiv(rhoPhi, Yi)
==
combustion->R(Yi)
+ parcels.SYi(i, Yi)
);

// solve
// fvScalarMatrix YiEqn
// (
// fvm::ddt(rhoYi[i])
// + fvc::div(rhoPhiYi)
// ==
// chemistry->RR(i)
// + parcels.SYi(i, rhoYi[i])
// fvm::ddt(rho, Yi)
// + mvConvection->fvmDiv(rhoPhi, Yi)
// ==
// combustion->R(Yi)
// + parcels.SYi(i, Yi)
// );

solve
(
fvm::ddt(rhoYi[i])
+ fvc::div(rhoPhiYi[i])
==
chemistry->RR(i)
+ parcels.SYi(i, rhoYi[i])
);

// Yi=rhoYi[i]/rho;
// Yi.max(0.0);
Yi=rhoYi[i]/rho;
Yi.max(0.0);

if (!inviscid)
{
const surfaceScalarField phiUc = linearInterpolate(sumYDiffError) & mesh.Sf();
tmp<volScalarField> DEff = chemistry->rhoD(i) + turbulence->mut()/Sct;

// original term
YiEqn -=
(
turbName == "laminar"
? (fvm::laplacian(DEff(), Yi) - mvConvection->fvmDiv(phiUc, Yi))
: fvm::laplacian(DEff(), Yi)
);
// YiEqn -=
// (
// turbName == "laminar"
// ? (fvm::laplacian(DEff(), Yi) - mvConvection->fvmDiv(phiUc, Yi))
// : fvm::laplacian(DEff(), Yi)
// );


// fvScalarMatrix YiEqn
// (
// fvm::ddt(rho, Yi) - fvc::ddt(rho, Yi)
// -
// (
// turbName == "laminar"
// ? (fvm::laplacian(DEff(), Yi) - mvConvection->fvmDiv(phiUc, Yi))
// : fvm::laplacian(DEff(), Yi)
// )
// );
fvScalarMatrix YiEqn
(
fvm::ddt(rho, Yi) - fvc::ddt(rho, Yi)
-
(
turbName == "laminar"
? (fvm::laplacian(DEff(), Yi) - mvConvection->fvmDiv(phiUc, Yi))
: fvm::laplacian(DEff(), Yi)
)
);

// YiEqn.relax();
YiEqn.relax();

// YiEqn.solve("Yi");
YiEqn.solve("Yi");

// Yi.max(0.0);
Yi.max(0.0);
}
// original term
YiEqn.relax();
YiEqn.solve("Yi");
Yi.max(0.0);
// YiEqn.relax();
// YiEqn.solve("Yi");
// Yi.max(0.0);
rhoYi[i] = rho*Yi;

}
Expand Down
10 changes: 10 additions & 0 deletions src/fluxSchemes/AUSMDV/AUSMDV.C
Original file line number Diff line number Diff line change
Expand Up @@ -73,13 +73,16 @@ void Foam::fluxSchemes::AUSMDV::createSavedFields()
void Foam::fluxSchemes::AUSMDV::calculateFluxes
(
const scalar& rhoOwn, const scalar& rhoNei,
const scalarList& rhoYiOwn,
const scalarList& rhoYiNei,
const vector& UOwn, const vector& UNei,
const scalar& eOwn, const scalar& eNei,
const scalar& pOwn, const scalar& pNei,
const scalar& cOwn, const scalar& cNei,
const vector& Sf,
scalar& phi,
scalar& rhoPhi,
scalarList& rhoPhiYi,
vector& rhoUPhi,
scalar& rhoEPhi,
const label facei, const label patchi
Expand Down Expand Up @@ -129,6 +132,13 @@ void Foam::fluxSchemes::AUSMDV::calculateFluxes
- (1 - caseA*caseB)*(caseA*0.125*(UvNei - cNei - UvOwn + cOwn)*(rhoNei - rhoOwn)*magSf
+ (1 - caseA)*caseB*0.125*(UvNei + cNei - UvOwn - cOwn)*(rhoNei - rhoOwn)*magSf);

forAll(rhoPhiYi,i)
{
rhoPhiYi[i] = (uPlus*rhoYiOwn[i] + uMinus*rhoYiNei[i])*magSf
- (1 - caseA*caseB)*(caseA*0.125*(UvNei - cNei - UvOwn + cOwn)*(rhoYiNei[i] - rhoYiOwn[i])*magSf
+ (1 - caseA)*caseB*0.125*(UvNei + cNei - UvOwn - cOwn)*(rhoYiNei[i] - rhoYiOwn[i])*magSf);
}

vector AUSMV((uPlus*rhoOwn*UOwn + uMinus*rhoNei*UNei)*magSf);
vector AUSMD(0.5*(rhoPhi*(UOwn+UNei) - mag(rhoPhi)*(UNei-UOwn)));

Expand Down
3 changes: 3 additions & 0 deletions src/fluxSchemes/AUSMDV/AUSMDV.H
Original file line number Diff line number Diff line change
Expand Up @@ -74,13 +74,16 @@ class AUSMDV
virtual void calculateFluxes
(
const scalar& rhoOwn, const scalar& rhoNei,
const scalarList& rhoYiOwn,
const scalarList& rhoYiNei,
const vector& UOwn, const vector& UNei,
const scalar& eOwn, const scalar& eNei,
const scalar& pOwn, const scalar& pNei,
const scalar& cOwn, const scalar& cNei,
const vector& Sf,
scalar& phi,
scalar& rhoPhi,
scalarList& rhoPhiYi,
vector& rhoUPhi,
scalar& rhoEPhi,
const label facei, const label patchi = -1
Expand Down
23 changes: 23 additions & 0 deletions src/fluxSchemes/HLLC/HLLC.C
Original file line number Diff line number Diff line change
Expand Up @@ -70,13 +70,16 @@ void Foam::fluxSchemes::HLLC::createSavedFields()
void Foam::fluxSchemes::HLLC::calculateFluxes
(
const scalar& rhoOwn, const scalar& rhoNei,
const scalarList& rhoYiOwn,
const scalarList& rhoYiNei,
const vector& UOwn, const vector& UNei,
const scalar& eOwn, const scalar& eNei,
const scalar& pOwn, const scalar& pNei,
const scalar& cOwn, const scalar& cNei,
const vector& Sf,
scalar& phi,
scalar& rhoPhi,
scalarList& rhoPhiYi,
vector& rhoUPhi,
scalar& rhoEPhi,
const label facei, const label patchi
Expand Down Expand Up @@ -142,6 +145,10 @@ void Foam::fluxSchemes::HLLC::calculateFluxes
// this->save(facei, patchi, UOwn, Uf_);
phi = UvOwn;
rhoPhi = rhoOwn*UvOwn;
forAll(rhoPhiYi,i)
{
rhoPhiYi[i] = rhoYiOwn[i]*UvOwn;
}
p = pOwn;
rhoUPhi = rhoUPhiOwn;
rhoEPhi = rhoEPhiOwn;
Expand All @@ -160,6 +167,10 @@ void Foam::fluxSchemes::HLLC::calculateFluxes
// );
phi = SStar;
rhoPhi = phi*rhoOwn*(SOwn - UvOwn)/dS;
forAll(rhoPhiYi,i)
{
rhoPhiYi[i] = phi*rhoYiOwn[i]*(SOwn - UvOwn)/dS;
}
p = 0.5*(pStarNei + pStarOwn);
rhoUPhi =
(SStar*(SOwn*rhoUOwn - rhoUPhiOwn) + SOwn*pStarOwn*normal)/dS;
Expand All @@ -179,6 +190,10 @@ void Foam::fluxSchemes::HLLC::calculateFluxes
// );
phi = SStar;
rhoPhi = phi*rhoNei*(SNei - UvNei)/dS;
forAll(rhoPhiYi,i)
{
rhoPhiYi[i] = phi*rhoYiNei[i]*(SNei - UvNei)/dS;
}
p = 0.5*(pStarNei + pStarOwn);
rhoUPhi =
(SStar*(SNei*rhoUNei - rhoUPhiNei) + SNei*pStarNei*normal)/dS;
Expand All @@ -189,12 +204,20 @@ void Foam::fluxSchemes::HLLC::calculateFluxes
// this->save(facei, patchi, UNei, Uf_);
phi = UvNei;
rhoPhi = rhoNei*UvNei;
forAll(rhoPhiYi,i)
{
rhoPhiYi[i] = rhoYiNei[i]*UvNei;
}
p = pNei;
rhoUPhi = rhoUPhiNei;
rhoEPhi = rhoEPhiNei;
}
phi *= magSf;
rhoPhi *= magSf;
forAll(rhoPhiYi,i)
{
rhoPhiYi[i] *= magSf;
}
rhoUPhi *= magSf;
rhoEPhi *= magSf;
rhoEPhi += vMesh*magSf*p;
Expand Down
3 changes: 3 additions & 0 deletions src/fluxSchemes/HLLC/HLLC.H
Original file line number Diff line number Diff line change
Expand Up @@ -81,13 +81,16 @@ class HLLC
virtual void calculateFluxes
(
const scalar& rhoOwn, const scalar& rhoNei,
const scalarList& rhoYiOwn,
const scalarList& rhoYiNei,
const vector& UOwn, const vector& UNei,
const scalar& eOwn, const scalar& eNei,
const scalar& pOwn, const scalar& pNei,
const scalar& cOwn, const scalar& cNei,
const vector& Sf,
scalar& phi,
scalar& rhoPhi,
scalarList& rhoPhiYi,
vector& rhoUPhi,
scalar& rhoEPhi,
const label facei, const label patchi = -1
Expand Down
Loading
Loading