diff --git a/src/thermophysicalModels/radiation/radiationModel/P1/P1.C b/src/thermophysicalModels/radiation/radiationModel/P1/P1.C index 5a7e5982789642aa1d02f85456fed0bf701adb12..dee96292927b8c0032fd5ab2400cfb77c1b78611 100644 --- a/src/thermophysicalModels/radiation/radiationModel/P1/P1.C +++ b/src/thermophysicalModels/radiation/radiationModel/P1/P1.C @@ -50,12 +50,9 @@ namespace Foam } } -// Radiation solver iteration -Foam::label Foam::radiation::P1::iterRadId = pTraits<label>::one; // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // -// Construct from components Foam::radiation::P1::P1(const volScalarField& T) : radiationModel(typeName, T), @@ -136,14 +133,8 @@ bool Foam::radiation::P1::read() } -void Foam::radiation::P1::correct() +void Foam::radiation::P1::calculate() { - if (!radiation_ || !(iterRadId == nFlowIterPerRadIter_)) - { - iterRadId++; - return; - } - a_ = absorptionEmission_->a(); e_ = absorptionEmission_->e(); E_ = absorptionEmission_->E(); @@ -171,8 +162,6 @@ void Foam::radiation::P1::correct() == - 4.0*(e_*radiation::sigmaSB*pow4(T_) + E_) ); - - iterRadId = pTraits<label>::one; } diff --git a/src/thermophysicalModels/radiation/radiationModel/P1/P1.H b/src/thermophysicalModels/radiation/radiationModel/P1/P1.H index b436a8509f5aff0b623b5c9445aa8e72a59a4bdd..c0249e28dbd592828a21fc78e18a2d4a501bd1eb 100644 --- a/src/thermophysicalModels/radiation/radiationModel/P1/P1.H +++ b/src/thermophysicalModels/radiation/radiationModel/P1/P1.H @@ -59,7 +59,6 @@ class P1 : public radiationModel { - // Private data //- Incident radiation / [W/m2] @@ -86,10 +85,6 @@ class P1 public: - // Static data members - - static label iterRadId; - //- Runtime type information TypeName("P1"); @@ -108,10 +103,10 @@ public: // Edit - //- Update radiationSource varible - void correct(); + //- Solve radiation equation(s) + void calculate(); - //- Read radiationProperties dictionary + //- Read radiation properties dictionary bool read(); diff --git a/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.C b/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.C index 4b57942cb7723456b2928f8dc2ea48bea907bc68..222da749ad7c7a35b4057a53b715e5ab8daaf8e5 100644 --- a/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.C +++ b/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.C @@ -50,9 +50,6 @@ namespace Foam } } -// Radiation solver iterator counter -Foam::label Foam::radiation::fvDOM::iterRadId = pTraits<label>::one; - // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // @@ -287,19 +284,13 @@ bool Foam::radiation::fvDOM::read() } -void Foam::radiation::fvDOM::correct() +void Foam::radiation::fvDOM::calculate() { - if (!radiation_ || !(iterRadId == nFlowIterPerRadIter_)) - { - iterRadId++; - return; - } - absorptionEmission_->correct(a_, aj_); updateBlackBodyEmission(); - scalar maxResidual = 0; + scalar maxResidual = 0.0; label radIter = 0; do { @@ -316,8 +307,6 @@ void Foam::radiation::fvDOM::correct() } while(maxResidual > convergence_); updateG(); - - iterRadId = pTraits<label>::one; } diff --git a/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.H b/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.H index 79899701778b5418ad96b7bbecae829859340f92..3c25cce890d25e3b486b3148bdebc1615f3a81b5 100644 --- a/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.H +++ b/src/thermophysicalModels/radiation/radiationModel/fvDOM/fvDOM/fvDOM.H @@ -136,11 +136,6 @@ class fvDOM public: - // Static data members - - static label iterRadId; - - //- Runtime type information TypeName("fvDOM"); @@ -159,10 +154,10 @@ public: // Edit - //- Update radiation source variable - void correct(); + //- Solve radiation equation(s) + void calculate(); - //- Read radiationProperties dictionary + //- Read radiation properties dictionary bool read(); //- Update G and calculate total heat flux on boundary diff --git a/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.C b/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.C index 6d387ce4e74ef05c013025d6f6a34041e4dfce29..062b07a85dee63a8cf9f4c385b93856a1b20d549 100644 --- a/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.C +++ b/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.C @@ -117,7 +117,7 @@ Foam::radiation::radiativeIntensityRay::radiativeIntensityRay { IOobject IHeader ( - "Ilambda_" + name(rayId) + "_" + name(i), + "ILambda_" + name(rayId) + "_" + name(i), mesh_.time().timeName(), mesh_, IOobject::MUST_READ, @@ -170,7 +170,7 @@ Foam::radiation::radiativeIntensityRay::~radiativeIntensityRay() Foam::scalar Foam::radiation::radiativeIntensityRay::correct() { // reset boundary heat flux to zero - Qr_ = dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0); + Qr_ = dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0); scalar maxResidual = -GREAT; diff --git a/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.H b/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.H index 0f26639064addb935697a51a9ac0bf92ba600730..27fbf1bc71ca8bc6023a997b31a701b201d4b1c1 100644 --- a/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.H +++ b/src/thermophysicalModels/radiation/radiationModel/fvDOM/radiativeIntensityRay/radiativeIntensityRay.H @@ -96,6 +96,9 @@ class radiativeIntensityRay //- List of pointers to radiative intensity fields for given wavelengths PtrList<volScalarField> ILambda_; + //- Global ray id - incremented in constructor + static label rayId; + // Private member functions @@ -108,11 +111,6 @@ class radiativeIntensityRay public: - // Static data members - - static label rayId; - - // Constructors //- Construct form components diff --git a/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.C b/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.C index ca814fe305a7fca33ce3fe564e35b49ffb12dfcd..fc6fbc549a3e3237ba7238436db5dbc8f93e2d8f 100644 --- a/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.C +++ b/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.C @@ -50,10 +50,9 @@ namespace Foam // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // -// Construct from components Foam::radiation::noRadiation::noRadiation(const volScalarField& T) : - radiationModel(typeName, T) + radiationModel(T) {} @@ -71,7 +70,7 @@ bool Foam::radiation::noRadiation::read() } -void Foam::radiation::noRadiation::correct() +void Foam::radiation::noRadiation::calculate() { // Do nothing } diff --git a/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.H b/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.H index 2c0e06dfde35bddb6a92492d345e4c7b1ebbdff4..2527ab1211fac521aa8d55b8fb4a4b90ca873316 100644 --- a/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.H +++ b/src/thermophysicalModels/radiation/radiationModel/noRadiation/noRadiation.H @@ -54,7 +54,6 @@ class noRadiation : public radiationModel { - // Private member functions //- Disallow default bitwise copy construct @@ -76,7 +75,7 @@ public: noRadiation(const volScalarField& T); - // Destructor + //- Destructor virtual ~noRadiation(); @@ -84,8 +83,8 @@ public: // Edit - //- Update radiation source - void correct(); + //- Solve radiation equation(s) + void calculate(); //- Read radiationProperties dictionary bool read(); diff --git a/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.C b/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.C index 28e10cd84eca9133bdc7e7c1b41c39fd4e62eebd..131221c308e3a74dcecf93b88af6a0bb45f9d872 100644 --- a/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.C +++ b/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.C @@ -43,6 +43,30 @@ namespace Foam // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // +Foam::radiation::radiationModel::radiationModel(const volScalarField& T) +: + IOdictionary + ( + IOobject + ( + "radiationProperties", + T.time().constant(), + T.mesh().objectRegistry::db(), + IOobject::MUST_READ, + IOobject::NO_WRITE + ) + ), + mesh_(T.mesh()), + time_(T.time()), + T_(T), + radiation_(false), + coeffs_(dictionary::null), + solverFreq_(0), + absorptionEmission_(NULL), + scatter_(NULL) +{} + + Foam::radiation::radiationModel::radiationModel ( const word& type, @@ -60,14 +84,17 @@ Foam::radiation::radiationModel::radiationModel IOobject::NO_WRITE ) ), - T_(T), mesh_(T.mesh()), + time_(T.time()), + T_(T), radiation_(lookup("radiation")), coeffs_(subDict(type + "Coeffs")), - nFlowIterPerRadIter_(readLabel(lookup("nFlowIterPerRadIter"))), + solverFreq_(readLabel(lookup("solverFreq"))), absorptionEmission_(absorptionEmissionModel::New(*this, mesh_)), scatter_(scatterModel::New(*this, mesh_)) -{} +{ + solverFreq_ = max(1, solverFreq_); +} // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * // @@ -94,6 +121,20 @@ bool Foam::radiation::radiationModel::read() } +void Foam::radiation::radiationModel::correct() +{ + if (!radiation_) + { + return; + } + + if ((time_.timeIndex() == 0) || (time_.timeIndex() % solverFreq_ == 0)) + { + calculate(); + } +} + + Foam::tmp<Foam::fvScalarMatrix> Foam::radiation::radiationModel::Sh ( basicThermo& thermo diff --git a/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.H b/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.H index 31e5dd25198345a74f7839cde4248a96b7cc4d28..4d01dee1cf158d11c2ce49b21cf3d16aa3934d6c 100644 --- a/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.H +++ b/src/thermophysicalModels/radiation/radiationModel/radiationModel/radiationModel.H @@ -73,21 +73,24 @@ protected: // Protected data + //- Reference to the mesh database + const fvMesh& mesh_; + + //- Reference to the time database + const Time& time_; + //- Reference to the temperature field const volScalarField& T_; - //- Reference to the mesh - const fvMesh& mesh_; - //- Model specific dictionary input parameters Switch radiation_; //- Radiation model dictionary dictionary coeffs_; - //- Number of iterations in the flow solver per radiation solver - // iteration - label nFlowIterPerRadIter_; + //- Radiation solver frequency - number flow solver iterations per + // radiation solver iteration + label solverFreq_; // References to the radiation sub-models @@ -132,6 +135,9 @@ public: // Constructors + //- Null constructor + radiationModel(const volScalarField& T); + //- Construct from components radiationModel(const word& type, const volScalarField& T); @@ -142,7 +148,7 @@ public: static autoPtr<radiationModel> New(const volScalarField& T); - // Destructor + //- Destructor virtual ~radiationModel(); @@ -150,8 +156,11 @@ public: // Edit + //- Main update/correction routine + virtual void correct(); + //- Solve radiation equation(s) - virtual void correct() = 0; + virtual void calculate() = 0; //- Read radiationProperties dictionary virtual bool read() = 0;