From e08f5bca5e35bea29beb2118826a846db6784795 Mon Sep 17 00:00:00 2001 From: andy <andy> Date: Mon, 15 Oct 2012 10:45:55 +0100 Subject: [PATCH] ENH: Added new contructor and selection table for radiation models --- .../radiationModel/radiationModel.C | 64 +++++++++++++++++++ .../radiationModel/radiationModel.H | 55 ++++++++++++---- .../radiationModel/radiationModelNew.C | 37 +++++++++-- 3 files changed, 140 insertions(+), 16 deletions(-) diff --git a/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.C b/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.C index 6a20fbcf5cd..57edaff7962 100644 --- a/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.C +++ b/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.C @@ -35,6 +35,7 @@ namespace Foam namespace radiation { defineTypeNameAndDebug(radiationModel, 0); + defineRunTimeSelectionTable(radiationModel, T); defineRunTimeSelectionTable(radiationModel, dictionary); } } @@ -67,6 +68,36 @@ Foam::radiation::radiationModel::radiationModel(const volScalarField& T) {} +Foam::radiation::radiationModel::radiationModel +( + const dictionary& dict, + const volScalarField& T +) +: + IOdictionary + ( + IOobject + ( + "radiationProperties", + T.time().constant(), + T.mesh(), + IOobject::NO_READ, + IOobject::NO_WRITE + ), + dict + ), + mesh_(T.mesh()), + time_(T.time()), + T_(T), + radiation_(false), + coeffs_(dictionary::null), + solverFreq_(0), + firstIter_(true), + absorptionEmission_(NULL), + scatter_(NULL) +{} + + Foam::radiation::radiationModel::radiationModel ( const word& type, @@ -98,6 +129,39 @@ Foam::radiation::radiationModel::radiationModel } +Foam::radiation::radiationModel::radiationModel +( + const word& type, + const dictionary& dict, + const volScalarField& T +) +: + IOdictionary + ( + IOobject + ( + "radiationProperties", + T.time().constant(), + T.mesh(), + IOobject::NO_READ, + IOobject::NO_WRITE + ), + dict + ), + mesh_(T.mesh()), + time_(T.time()), + T_(T), + radiation_(lookup("radiation")), + coeffs_(subDict(type + "Coeffs")), + solverFreq_(lookupOrDefault<label>("solverFreq", 1)), + firstIter_(true), + absorptionEmission_(absorptionEmissionModel::New(*this, mesh_)), + scatter_(scatterModel::New(*this, mesh_)) +{ + solverFreq_ = max(1, solverFreq_); +} + + // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * // Foam::radiation::radiationModel::~radiationModel() diff --git a/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.H b/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.H index e8aa9f5778e..d9c48053e71 100644 --- a/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.H +++ b/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModel.H @@ -35,6 +35,7 @@ Description SourceFiles radiationModel.C + radiationModelNew.C \*---------------------------------------------------------------------------*/ @@ -123,16 +124,28 @@ public: // Declare runtime constructor selection table - declareRunTimeSelectionTable - ( - autoPtr, - radiationModel, - dictionary, - ( - const volScalarField& T - ), - (T) - ); + declareRunTimeSelectionTable + ( + autoPtr, + radiationModel, + T, + ( + const volScalarField& T + ), + (T) + ); + + declareRunTimeSelectionTable + ( + autoPtr, + radiationModel, + dictionary, + ( + const dictionary& dict, + const volScalarField& T + ), + (dict, T) + ); // Constructors @@ -140,14 +153,32 @@ public: //- Null constructor radiationModel(const volScalarField& T); + //- Construct with dictionary + radiationModel(const dictionary& dict, const volScalarField& T); + //- Construct from components radiationModel(const word& type, const volScalarField& T); + //- Construct from components + radiationModel + ( + const word& type, + const dictionary& dict, + const volScalarField& T + ); + // Selectors - //- Return a reference to the selected radiation model - static autoPtr<radiationModel> New(const volScalarField& T); + //- Return a reference to the selected radiation model + static autoPtr<radiationModel> New(const volScalarField& T); + + //- Return a reference to the selected radiation model + static autoPtr<radiationModel> New + ( + const dictionary& dict, + const volScalarField& T + ); //- Destructor diff --git a/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModelNew.C b/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModelNew.C index 8161172f25e..ae13f11d1a5 100644 --- a/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModelNew.C +++ b/src/thermophysicalModels/radiationModels/radiationModel/radiationModel/radiationModelNew.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -23,7 +23,6 @@ License \*---------------------------------------------------------------------------*/ -#include "error.H" #include "radiationModel.H" // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // @@ -53,6 +52,36 @@ Foam::radiation::radiationModel::New Info<< "Selecting radiationModel " << modelType << endl; + TConstructorTable::iterator cstrIter = + TConstructorTablePtr_->find(modelType); + + if (cstrIter == TConstructorTablePtr_->end()) + { + FatalErrorIn + ( + "radiationModel::New(const volScalarField&)" + ) << "Unknown radiationModel type " + << modelType << nl << nl + << "Valid radiationModel types are:" << nl + << TConstructorTablePtr_->sortedToc() + << exit(FatalError); + } + + return autoPtr<radiationModel>(cstrIter()(T)); +} + + +Foam::autoPtr<Foam::radiation::radiationModel> +Foam::radiation::radiationModel::New +( + const dictionary& dict, + const volScalarField& T +) +{ + const word modelType(dict.lookup("radiationModel")); + + Info<< "Selecting radiationModel " << modelType << endl; + dictionaryConstructorTable::iterator cstrIter = dictionaryConstructorTablePtr_->find(modelType); @@ -60,7 +89,7 @@ Foam::radiation::radiationModel::New { FatalErrorIn ( - "radiationModel::New(const volScalarField&)" + "radiationModel::New(const dictionary&, const volScalarField&)" ) << "Unknown radiationModel type " << modelType << nl << nl << "Valid radiationModel types are:" << nl @@ -68,7 +97,7 @@ Foam::radiation::radiationModel::New << exit(FatalError); } - return autoPtr<radiationModel>(cstrIter()(T)); + return autoPtr<radiationModel>(cstrIter()(dict, T)); } -- GitLab