From 45cbffa4566860cb5127ed5dd2a52aa733ee4e32 Mon Sep 17 00:00:00 2001 From: andy <a.heather@opencfd.co.uk> Date: Mon, 25 Oct 2010 15:58:31 +0100 Subject: [PATCH] ENH: Updated cloud source treatments --- .../KinematicCloud/KinematicCloudI.H | 4 +- .../Templates/ReactingCloud/ReactingCloud.H | 10 ++ .../Templates/ReactingCloud/ReactingCloudI.H | 132 ++++++++++++++++++ .../Templates/ThermoCloud/ThermoCloud.C | 3 +- .../Templates/ThermoCloud/ThermoCloudI.H | 8 +- 5 files changed, 147 insertions(+), 10 deletions(-) diff --git a/src/lagrangian/intermediate/clouds/Templates/KinematicCloud/KinematicCloudI.H b/src/lagrangian/intermediate/clouds/Templates/KinematicCloud/KinematicCloudI.H index 94b816c8e7d..3643f52cd50 100644 --- a/src/lagrangian/intermediate/clouds/Templates/KinematicCloud/KinematicCloudI.H +++ b/src/lagrangian/intermediate/clouds/Templates/KinematicCloud/KinematicCloudI.H @@ -439,11 +439,11 @@ Foam::KinematicCloud<ParcelType>::SU(volVectorField& U) const if (solution_.sourceActive()) { - if (solution_.semiImplicit("UTrans")) + if (solution_.semiImplicit("U")) { return UTrans()/(mesh_.V()*this->db().time().deltaT()) - + fvm::SuSp(-UCoeff().dimensionedInternalField()/mesh_.V(), U) + + fvm::SuSp(-UCoeff()/mesh_.V(), U) + UCoeff()/mesh_.V()*U; } else diff --git a/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.H b/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.H index 539c787c102..c968a723f10 100644 --- a/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.H +++ b/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.H @@ -219,6 +219,13 @@ public: inline PtrList<DimensionedField<scalar, volMesh> >& rhoTrans(); + //- Return mass source term for specie i - specie eqn + inline tmp<fvScalarMatrix> SYi + ( + const label i, + volScalarField& Yi + ) const; + //- Return tmp mass source for field i - fully explicit inline tmp<DimensionedField<scalar, volMesh> > Srho(const label i) const; @@ -227,6 +234,9 @@ public: // - fully explicit inline tmp<DimensionedField<scalar, volMesh> > Srho() const; + //- Return total mass source term [kg/m3/s] + inline tmp<fvScalarMatrix> Srho(volScalarField& rho) const; + // Check diff --git a/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloudI.H b/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloudI.H index a94de260c20..b9c23839ec3 100644 --- a/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloudI.H +++ b/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloudI.H @@ -65,6 +65,62 @@ Foam::ReactingCloud<ParcelType>::rhoTrans() } +template<class ParcelType> +inline Foam::tmp<Foam::fvScalarMatrix> +Foam::ReactingCloud<ParcelType>::SYi +( + const label i, + volScalarField& Yi +) const +{ + if (this->solution().sourceActive()) + { + if (this->solution().semiImplicit("Yi")) + { + tmp<volScalarField> trhoTrans + ( + new volScalarField + ( + IOobject + ( + this->name() + "rhoTrans", + this->db().time().timeName(), + this->db(), + IOobject::NO_READ, + IOobject::NO_WRITE, + false + ), + this->mesh(), + dimensionedScalar("zero", dimMass/dimTime/dimVolume, 0.0) + ) + ); + + volScalarField& sourceField = trhoTrans(); + + sourceField.internalField() = + rhoTrans_[i]/(this->db().time().deltaT()*this->mesh().V()); + + const dimensionedScalar YiSMALL("YiSMALL", dimless, SMALL); + + return + fvm::Sp(neg(sourceField)*sourceField/(Yi + YiSMALL), Yi) + + pos(sourceField)*sourceField; + } + else + { + tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(Yi, dimMass/dimTime)); + fvScalarMatrix& fvm = tfvm(); + + fvm.source() = -rhoTrans_[i]/this->db().time().deltaT(); + + return tfvm; + } + } + + return tmp<fvScalarMatrix>(new fvScalarMatrix(Yi, dimMass/dimTime)); +} + + template<class ParcelType> inline Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::volMesh> > Foam::ReactingCloud<ParcelType>::Srho(const label i) const @@ -166,4 +222,80 @@ Foam::ReactingCloud<ParcelType>::Srho() const } +template<class ParcelType> +inline Foam::tmp<Foam::fvScalarMatrix> +Foam::ReactingCloud<ParcelType>::Srho(volScalarField& rho) const +{ + if (this->solution().sourceActive()) + { + if (this->solution().semiImplicit("rho")) + { + tmp<volScalarField> trhoTrans + ( + new volScalarField + ( + IOobject + ( + this->name() + "rhoTrans", + this->db().time().timeName(), + this->db(), + IOobject::NO_READ, + IOobject::NO_WRITE, + false + ), + this->mesh(), + dimensionedScalar("zero", dimMass/dimTime/dimVolume, 0.0) + ) + ); + + scalarField& sourceField = trhoTrans(); + + forAll(rhoTrans_, i) + { + sourceField += rhoTrans_[i]; + } + sourceField /= this->db().time().deltaT()*this->mesh().V(); + + return fvm::SuSp(trhoTrans()/rho, rho); + } + else + { + tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(rho, dimMass/dimTime)); + fvScalarMatrix& fvm = tfvm(); + + tmp<volScalarField> trhoTrans + ( + new volScalarField + ( + IOobject + ( + this->name() + "rhoTrans", + this->db().time().timeName(), + this->db(), + IOobject::NO_READ, + IOobject::NO_WRITE, + false + ), + this->mesh(), + dimensionedScalar("zero", dimMass, 0.0) + ) + ); + + scalarField& sourceField = trhoTrans(); + + forAll(rhoTrans_, i) + { + sourceField += rhoTrans_[i]; + } + + fvm.source() = -trhoTrans()/this->db().time().deltaT(); + + return tfvm; + } + } + + return tmp<fvScalarMatrix>(new fvScalarMatrix(rho, dimMass/dimTime)); +} + + // ************************************************************************* // diff --git a/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C b/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C index 69cde7fe7c0..61a19d47626 100644 --- a/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C +++ b/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C @@ -117,8 +117,7 @@ Foam::ThermoCloud<ParcelType>::ThermoCloud IOobject::AUTO_WRITE ), this->mesh(), - dimensionedScalar("zero", dimEnergy/dimTime/dimTemperature, 0.0), - zeroGradientFvPatchScalarField::typeName + dimensionedScalar("zero", dimEnergy/dimTime/dimTemperature, 0.0) ) ) diff --git a/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloudI.H b/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloudI.H index e1c6e88f14c..79e26fa1661 100644 --- a/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloudI.H +++ b/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloudI.H @@ -127,17 +127,13 @@ Foam::ThermoCloud<ParcelType>::Sh(volScalarField& hs) const if (this->solution().sourceActive()) { - if (this->solution().semiImplicit("hsTrans")) + if (this->solution().semiImplicit("hs")) { const volScalarField Cp = thermo_.thermo().Cp(); return hsTrans()/(this->mesh().V()*this->db().time().deltaT()) - + fvm::SuSp - ( - -hsCoeff().dimensionedInternalField()/(Cp*this->mesh().V()), - hs - ) + + fvm::SuSp(-hsCoeff()/(Cp*this->mesh().V()), hs) + hsCoeff()/(Cp*this->mesh().V())*hs; } else -- GitLab