From 3e5be0f42e89f11c7b505df1b3823a9a51f2e362 Mon Sep 17 00:00:00 2001 From: Thomas Lichtenegger Date: Mon, 26 Apr 2021 16:09:51 +0200 Subject: [PATCH] Added otherForceModels to PISO solver. --- applications/solvers/cfdemSolverPiso/UEqn.H | 3 +++ .../solvers/cfdemSolverPiso/createFields.H | 15 +++++++++++++++ .../solvers/cfdemSolverPisoScalar/createFields.H | 15 +++++++++++++++ 3 files changed, 33 insertions(+) diff --git a/applications/solvers/cfdemSolverPiso/UEqn.H b/applications/solvers/cfdemSolverPiso/UEqn.H index 2f26ac91..21165a36 100644 --- a/applications/solvers/cfdemSolverPiso/UEqn.H +++ b/applications/solvers/cfdemSolverPiso/UEqn.H @@ -1,8 +1,11 @@ +particleCloud.otherForces(fOther); + fvVectorMatrix UEqn ( fvm::ddt(voidfraction,U) - fvm::Sp(fvc::ddt(voidfraction),U) + fvm::div(phi,U) - fvm::Sp(fvc::div(phi),U) + particleCloud.divVoidfractionTau(U, voidfraction) + - fOther/rho == fvOptions(U) - fvm::Sp(Ksl/rho,U) diff --git a/applications/solvers/cfdemSolverPiso/createFields.H b/applications/solvers/cfdemSolverPiso/createFields.H index 36246301..72f12354 100644 --- a/applications/solvers/cfdemSolverPiso/createFields.H +++ b/applications/solvers/cfdemSolverPiso/createFields.H @@ -46,6 +46,21 @@ //dimensionedScalar("0", dimensionSet(1, -3, -1, 0, 0), 1.0) ); + Info<< "\nCreating body force field\n" << endl; + volVectorField fOther + ( + IOobject + ( + "fOther", + runTime.timeName(), + mesh, + IOobject::NO_READ, + IOobject::NO_WRITE + ), + mesh, + dimensionedVector("zero", dimensionSet(1,-2,-2,0,0,0,0), vector::zero) + ); + Info<< "\nReading voidfraction field voidfraction = (Vgas/Vparticle)\n" << endl; volScalarField voidfraction ( diff --git a/applications/solvers/cfdemSolverPisoScalar/createFields.H b/applications/solvers/cfdemSolverPisoScalar/createFields.H index f7779bec..9cd6599f 100644 --- a/applications/solvers/cfdemSolverPisoScalar/createFields.H +++ b/applications/solvers/cfdemSolverPisoScalar/createFields.H @@ -46,6 +46,21 @@ //dimensionedScalar("0", dimensionSet(0, 0, -1, 0, 0), 1.0) ); + Info<< "\nCreating body force field\n" << endl; + volVectorField fOther + ( + IOobject + ( + "fOther", + runTime.timeName(), + mesh, + IOobject::NO_READ, + IOobject::NO_WRITE + ), + mesh, + dimensionedVector("zero", dimensionSet(1,-2,-2,0,0,0,0), vector::zero) + ); + Info<< "\nReading voidfraction field voidfraction = (Vgas/Vparticle)\n" << endl; volScalarField voidfraction (