radiationModelNew.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "radiationModel.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const volScalarField& T
34 )
35 {
36  IOobject radIO
37  (
38  "radiationProperties",
39  T.time().constant(),
40  T.mesh(),
41  IOobject::MUST_READ_IF_MODIFIED,
42  IOobject::NO_WRITE,
43  false
44  );
45 
46  word modelType("none");
47  if (radIO.typeHeaderOk<IOdictionary>(false))
48  {
49  IOdictionary(radIO).lookup("radiationModel") >> modelType;
50  }
51  else
52  {
53  Info<< "Radiation model not active: radiationProperties not found"
54  << endl;
55  }
56 
57  Info<< "Selecting radiationModel " << modelType << endl;
58 
59  TConstructorTable::iterator cstrIter =
60  TConstructorTablePtr_->find(modelType);
61 
62  if (cstrIter == TConstructorTablePtr_->end())
63  {
65  << "Unknown radiationModel type "
66  << modelType << nl << nl
67  << "Valid radiationModel types are:" << nl
68  << TConstructorTablePtr_->sortedToc()
69  << exit(FatalError);
70  }
71 
72  return autoPtr<radiationModel>(cstrIter()(T));
73 }
74 
75 
77 (
78  const dictionary& dict,
79  const volScalarField& T
80 )
81 {
82  const word modelType(dict.lookup("radiationModel"));
83 
84  Info<< "Selecting radiationModel " << modelType << endl;
85 
86  dictionaryConstructorTable::iterator cstrIter =
87  dictionaryConstructorTablePtr_->find(modelType);
88 
89  if (cstrIter == dictionaryConstructorTablePtr_->end())
90  {
92  << "Unknown radiationModel type "
93  << modelType << nl << nl
94  << "Valid radiationModel types are:" << nl
95  << dictionaryConstructorTablePtr_->sortedToc()
96  << exit(FatalError);
97  }
98 
99  return autoPtr<radiationModel>(cstrIter()(dict, T));
100 }
101 
102 
103 // ************************************************************************* //
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
static const char nl
Definition: Ostream.H:265
const Mesh & mesh() const
Return mesh.
const volScalarField & T
const Time & time() const
Return time.
Definition: IOobject.C:360
messageStream Info
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583