Skip to content
Snippets Groups Projects
Commit c9e57ad9 authored by Andrew Heather's avatar Andrew Heather
Browse files

updated to reflect change in lagrangian phase source term names

parent 30c7e8aa
Branches
Tags
No related merge requests found
...@@ -4,7 +4,6 @@ EXE_INC = \ ...@@ -4,7 +4,6 @@ EXE_INC = \
-I$(LIB_SRC)/turbulenceModels/compressible/turbulenceModel \ -I$(LIB_SRC)/turbulenceModels/compressible/turbulenceModel \
-I$(LIB_SRC)/lagrangian/basic/lnInclude \ -I$(LIB_SRC)/lagrangian/basic/lnInclude \
-I$(LIB_SRC)/lagrangian/intermediate/lnInclude \ -I$(LIB_SRC)/lagrangian/intermediate/lnInclude \
-I$(LIB_SRC)/lagrangian/coalCombustion/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/pdfs/lnInclude \ -I$(LIB_SRC)/thermophysicalModels/pdfs/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/specie/lnInclude \ -I$(LIB_SRC)/thermophysicalModels/specie/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \ -I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
...@@ -19,7 +18,6 @@ EXE_INC = \ ...@@ -19,7 +18,6 @@ EXE_INC = \
-I$(LIB_SRC)/ODE/lnInclude -I$(LIB_SRC)/ODE/lnInclude
EXE_LIBS = \ EXE_LIBS = \
-L$(FOAM_USER_LIBBIN) \
-lfiniteVolume \ -lfiniteVolume \
-lmeshTools \ -lmeshTools \
-lcompressibleRASModels \ -lcompressibleRASModels \
......
...@@ -5,7 +5,7 @@ ...@@ -5,7 +5,7 @@
+ turbulence->divDevRhoReff(U) + turbulence->divDevRhoReff(U)
== ==
rho.dimensionedInternalField()*g rho.dimensionedInternalField()*g
+ reactingParcels.SU1() + parcels.SU()
); );
UEqn.relax(); UEqn.relax();
......
tmp<fv::convectionScheme<scalar> > mvConvection tmp<fv::convectionScheme<scalar> > mvConvection
( (
fv::convectionScheme<scalar>::New fv::convectionScheme<scalar>::New
...@@ -26,7 +25,7 @@ tmp<fv::convectionScheme<scalar> > mvConvection ...@@ -26,7 +25,7 @@ tmp<fv::convectionScheme<scalar> > mvConvection
+ mvConvection->fvmDiv(phi, Yi) + mvConvection->fvmDiv(phi, Yi)
- fvm::laplacian(turbulence->muEff(), Yi) - fvm::laplacian(turbulence->muEff(), Yi)
== ==
reactingParcels.Srho1(i) parcels.Srho(i)
+ kappa*chemistry.RR(i)().dimensionedInternalField() + kappa*chemistry.RR(i)().dimensionedInternalField()
); );
...@@ -41,5 +40,4 @@ tmp<fv::convectionScheme<scalar> > mvConvection ...@@ -41,5 +40,4 @@ tmp<fv::convectionScheme<scalar> > mvConvection
Y[inertIndex] = scalar(1) - Yt; Y[inertIndex] = scalar(1) - Yt;
Y[inertIndex].max(0.0); Y[inertIndex].max(0.0);
} }
Info<< "\nConstructing interpolation" << endl;
Info << "\nConstructing gas properties" << endl; Info << "\nConstructing gas properties" << endl;
/* /*
PtrList<specieConstProperties> gasProperties(Y.size()); PtrList<specieConstProperties> gasProperties(Y.size());
...@@ -31,7 +29,7 @@ forAll(gasProperties, i) ...@@ -31,7 +29,7 @@ forAll(gasProperties, i)
} }
Info<< "\nConstructing reacting cloud" << endl; Info<< "\nConstructing reacting cloud" << endl;
basicReactingCloud reactingParcels basicReactingCloud parcels
( (
"reactingCloud1", "reactingCloud1",
rho, rho,
......
...@@ -28,51 +28,6 @@ ...@@ -28,51 +28,6 @@
thermo->rho() thermo->rho()
); );
// lagrangian coal density field
/* volScalarField rholc
(
IOobject
(
"rholc",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("zero", dimensionSet(1, -3, 0, 0, 0, 0, 0), 0.0)
);
// lagrangian limestone density field
volScalarField rhols
(
IOobject
(
"rhols",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("zero", dimensionSet(1, -3, 0, 0, 0, 0, 0), 0.0)
);
// lagrangian total density field
volScalarField rhol
(
IOobject
(
"rhol",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("zero", dimensionSet(1, -3, 0, 0, 0, 0, 0), 0.0)
);*/
Info<< "\nReading field U\n" << endl; Info<< "\nReading field U\n" << endl;
volVectorField U volVectorField U
( (
...@@ -133,9 +88,3 @@ ...@@ -133,9 +88,3 @@
fields.add(Y[i]); fields.add(Y[i]);
} }
fields.add(h); fields.add(h);
Info<< "Creating radiation model\n" << endl;
autoPtr<radiation::radiationModel> radiation
(
radiation::radiationModel::New(T)
);
...@@ -6,7 +6,7 @@ ...@@ -6,7 +6,7 @@
- fvm::laplacian(turbulence->alphaEff(), h) - fvm::laplacian(turbulence->alphaEff(), h)
== ==
DpDt DpDt
+ reactingParcels.Sh1() + parcels.Sh()
+ radiation->Sh(thermo()) + radiation->Sh(thermo())
); );
......
...@@ -23,7 +23,7 @@ if (transonic) ...@@ -23,7 +23,7 @@ if (transonic)
+ fvm::div(phid, p) + fvm::div(phid, p)
- fvm::laplacian(rho*rUA, p) - fvm::laplacian(rho*rUA, p)
== ==
reactingParcels.Srho1() parcels.Srho()
); );
pEqn.solve(); pEqn.solve();
...@@ -51,7 +51,7 @@ else ...@@ -51,7 +51,7 @@ else
+ fvc::div(phi) + fvc::div(phi)
- fvm::laplacian(rho*rUA, p) - fvm::laplacian(rho*rUA, p)
== ==
reactingParcels.Srho1() parcels.Srho()
); );
pEqn.solve(); pEqn.solve();
......
...@@ -23,8 +23,11 @@ License ...@@ -23,8 +23,11 @@ License
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Application Application
reactingParcelFoam
Description Description
Transient PISO solver for compressible, laminar or turbulent flow with
reacting Lagrangian parcels.
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
...@@ -42,55 +45,54 @@ Description ...@@ -42,55 +45,54 @@ Description
int main(int argc, char *argv[]) int main(int argc, char *argv[])
{ {
# include "setRootCase.H" #include "setRootCase.H"
# include "createTime.H" #include "createTime.H"
# include "createMesh.H" #include "createMesh.H"
# include "readChemistryProperties.H" #include "readChemistryProperties.H"
# include "readEnvironmentalProperties.H" #include "readEnvironmentalProperties.H"
# include "createFields.H" #include "createFields.H"
# include "createClouds.H" #include "createClouds.H"
# include "readPISOControls.H" #include "createRadiationModel.H"
# include "initContinuityErrs.H" #include "readPISOControls.H"
# include "readTimeControls.H" #include "initContinuityErrs.H"
# include "compressibleCourantNo.H" #include "readTimeControls.H"
# include "setInitialDeltaT.H" #include "compressibleCourantNo.H"
#include "setInitialDeltaT.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl; Info<< "\nStarting time loop\n" << endl;
while (runTime.run()) while (runTime.run())
{ {
# include "readTimeControls.H" #include "readTimeControls.H"
# include "readPISOControls.H" #include "readPISOControls.H"
# include "compressibleCourantNo.H" #include "compressibleCourantNo.H"
# include "setDeltaT.H" #include "setDeltaT.H"
runTime++; runTime++;
Info<< "Time = " << runTime.timeName() << nl << endl; Info<< "Time = " << runTime.timeName() << nl << endl;
Info << "Evolving reacting cloud" << endl; parcels.evolve();
reactingParcels.evolve();
reactingParcels.info(); parcels.info();
# include "chemistry.H" #include "chemistry.H"
# include "rhoEqn.H" #include "rhoEqn.H"
// --- PIMPLE loop // --- PIMPLE loop
for (int ocorr=1; ocorr<=nOuterCorr; ocorr++) for (int ocorr=1; ocorr<=nOuterCorr; ocorr++)
{ {
# include "UEqn.H" #include "UEqn.H"
# include "YEqn.H" #include "YEqn.H"
// --- PISO loop // --- PISO loop
for (int corr=1; corr<=nCorr; corr++) for (int corr=1; corr<=nCorr; corr++)
{ {
# include "hEqn.H" #include "hEqn.H"
# include "pEqn.H" #include "pEqn.H"
} }
Info<< "T gas min/max = " << min(T).value() << ", " Info<< "T gas min/max = " << min(T).value() << ", "
...@@ -103,7 +105,7 @@ int main(int argc, char *argv[]) ...@@ -103,7 +105,7 @@ int main(int argc, char *argv[])
if (runTime.write()) if (runTime.write())
{ {
# include "additionalOutput.H" #include "additionalOutput.H"
} }
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
......
...@@ -36,7 +36,7 @@ Description ...@@ -36,7 +36,7 @@ Description
fvm::ddt(rho) fvm::ddt(rho)
+ fvc::div(phi) + fvc::div(phi)
== ==
reactingParcels.Srho1() parcels.Srho()
); );
} }
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment