noThermo.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-2018 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 "noThermo.H"
28 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace thermalBaffleModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
42 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
50  return regionModel1D::read();
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 noThermo::noThermo(const word& modelType, const fvMesh& mesh)
57 :
58  thermalBaffleModel(mesh)
59 {}
60 
61 
62 noThermo::noThermo
63 (
64  const word& modelType,
65  const fvMesh& mesh,
66  const dictionary& dict
67 )
68 :
69  thermalBaffleModel(modelType, mesh, dict)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
76 {}
77 
78 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
79 
81 {}
82 
83 
85 {}
86 
87 
89 {
91  << "Cp field not available for " << type()
92  << abort(FatalError);
93 
94  return tmp<volScalarField>
95  (
96  new volScalarField
97  (
98  IOobject
99  (
100  "noThermo::Cp",
101  time().timeName(),
102  primaryMesh(),
105  false
106  ),
107  primaryMesh(),
109  )
110  );
111 }
112 
114 {
116  << "kappa field not available for " << type()
117  << abort(FatalError);
118  return volScalarField::null();
119 }
120 
121 
123 {
125  << "rho field not available for " << type()
126  << abort(FatalError);
127  return volScalarField::null();
128 }
129 
130 
132 {
134  << "K field not available for " << type()
135  << abort(FatalError);
136  return volScalarField::null();
137 }
138 
139 
141 {
143  << "T field not available for " << type()
144  << abort(FatalError);
145  return volScalarField::null();
146 }
147 
148 
150 {
152  << "T field not available for " << type()
153  << abort(FatalError);
154  return NullObjectRef<solidThermo>();
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace thermalBaffleModels
161 } // End namespace regionModels
162 } // End namespace Foam
163 
164 // ************************************************************************* //
dictionary dict
virtual void preEvolveRegion()
Pre-evolve film.
Definition: noThermo.C:80
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual const volScalarField & kappa() const
Return thermal conductivity [W/m/K].
Definition: noThermo.C:131
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
dynamicFvMesh & mesh
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
Definition: noThermo.C:149
addToRunTimeSelectionTable(thermalBaffleModel, noThermo, mesh)
A class for handling words, derived from string.
Definition: word.H:59
virtual const volScalarField & rho() const
Return density [Kg/m3].
Definition: noThermo.C:122
word timeName
Definition: getTimeIndex.H:3
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:48
virtual void evolveRegion()
Evolve the film equations.
Definition: noThermo.C:84
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
const dimensionSet dimEnergy
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
virtual const tmp< volScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: noThermo.C:88
Dummy surface pyrolysis model for &#39;none&#39;.
Definition: noThermo.H:53
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Read control parameters from dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual bool read()
Read control parameters from dictionary.
Definition: noThermo.C:48
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: noThermo.C:140
Namespace for OpenFOAM.
virtual const volScalarField & kappaRad() const
Return solid absortivity [1/m].
Definition: noThermo.C:113