heSolidThermo.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "heSolidThermo.H"
30#include "volFields.H"
31
32// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33
34template<class BasicSolidThermo, class MixtureType>
36{
37 scalarField& TCells = this->T_.primitiveFieldRef();
38
39 const scalarField& hCells = this->he_;
40 const scalarField& pCells = this->p_;
41 scalarField& rhoCells = this->rho_.primitiveFieldRef();
42 scalarField& alphaCells = this->alpha_.primitiveFieldRef();
43
44 forAll(TCells, celli)
45 {
46 const typename MixtureType::thermoType& mixture_ =
47 this->cellMixture(celli);
48
49 const typename MixtureType::thermoType& volMixture_ =
50 this->cellVolMixture(pCells[celli], TCells[celli], celli);
51
52 if (this->updateT())
53 {
54 TCells[celli] = mixture_.THE
55 (
56 hCells[celli],
57 pCells[celli],
58 TCells[celli]
59 );
60 }
61
62 rhoCells[celli] = volMixture_.rho(pCells[celli], TCells[celli]);
63
64 alphaCells[celli] =
65 volMixture_.kappa(pCells[celli], TCells[celli])
66 /
67 mixture_.Cpv(pCells[celli], TCells[celli]);
68 }
69
70 volScalarField::Boundary& pBf = this->p_.boundaryFieldRef();
71 volScalarField::Boundary& TBf = this->T_.boundaryFieldRef();
72 volScalarField::Boundary& rhoBf = this->rho_.boundaryFieldRef();
73 volScalarField::Boundary& heBf = this->he().boundaryFieldRef();
74 volScalarField::Boundary& alphaBf = this->alpha_.boundaryFieldRef();
75
76 forAll(this->T_.boundaryField(), patchi)
77 {
78 fvPatchScalarField& pp = pBf[patchi];
79 fvPatchScalarField& pT = TBf[patchi];
80 fvPatchScalarField& prho = rhoBf[patchi];
81 fvPatchScalarField& phe = heBf[patchi];
82 fvPatchScalarField& palpha = alphaBf[patchi];
83
84 if (pT.fixesValue())
85 {
86 forAll(pT, facei)
87 {
88 const typename MixtureType::thermoType& mixture_ =
89 this->patchFaceMixture(patchi, facei);
90
91 const typename MixtureType::thermoType& volMixture_ =
92 this->patchFaceVolMixture
93 (
94 pp[facei],
95 pT[facei],
96 patchi,
97 facei
98 );
99
100 phe[facei] = mixture_.HE(pp[facei], pT[facei]);
101 prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
102
103 palpha[facei] =
104 volMixture_.kappa(pp[facei], pT[facei])
105 / mixture_.Cpv(pp[facei], pT[facei]);
106 }
107 }
108 else
109 {
110 forAll(pT, facei)
111 {
112 const typename MixtureType::thermoType& mixture_ =
113 this->patchFaceMixture(patchi, facei);
114
115 const typename MixtureType::thermoType& volMixture_ =
116 this->patchFaceVolMixture
117 (
118 pp[facei],
119 pT[facei],
120 patchi,
121 facei
122 );
123
124 if (this->updateT())
125 {
126 pT[facei] = mixture_.THE(phe[facei], pp[facei] ,pT[facei]);
127 }
128
129 prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
130
131 palpha[facei] =
132 volMixture_.kappa(pp[facei], pT[facei])
133 / mixture_.Cpv(pp[facei], pT[facei]);
134 }
135 }
136 }
137
138 this->alpha_.correctBoundaryConditions();
139}
140
141// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142
143template<class BasicSolidThermo, class MixtureType>
146(
147 const fvMesh& mesh,
148 const word& phaseName
149)
150:
151 heThermo<BasicSolidThermo, MixtureType>(mesh, phaseName)
152{
153 calculate();
154 this->mu_ == dimensionedScalar(this->mu_.dimensions(), Zero);
155 this->psi_ == dimensionedScalar(this->psi_.dimensions(), Zero);
156}
157
158
159template<class BasicSolidThermo, class MixtureType>
162(
163 const fvMesh& mesh,
164 const dictionary& dict,
165 const word& phaseName
166)
167:
168 heThermo<BasicSolidThermo, MixtureType>(mesh, dict, phaseName)
169{
170 calculate();
171 this->mu_ == dimensionedScalar(this->mu_.dimensions(), Zero);
172 this->psi_ == dimensionedScalar(this->psi_.dimensions(), Zero);
173}
174
175
176template<class BasicSolidThermo, class MixtureType>
179(
180 const fvMesh& mesh,
181 const word& phaseName,
182 const word& dictName
183)
184:
185 heThermo<BasicSolidThermo, MixtureType>(mesh, phaseName, dictName)
186{
187 calculate();
188
189 // TBD. initialise psi, mu (at heThermo level) since these do not
190 // get initialised. Move to heThermo constructor?
191 this->mu_ == dimensionedScalar(this->mu_.dimensions(), Zero);
192 this->psi_ == dimensionedScalar(this->psi_.dimensions(), Zero);
193}
194
195
196// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
197
198template<class BasicSolidThermo, class MixtureType>
200{}
201
202
203// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204
205template<class BasicSolidThermo, class MixtureType>
207{
209
210 calculate();
211
212 DebugInfo << " Finished" << endl;
213}
214
215
216template<class BasicSolidThermo, class MixtureType>
219{
220 const fvMesh& mesh = this->T_.mesh();
221
223 (
225 (
227 (
228 "Kappa",
229 mesh.time().timeName(),
230 mesh,
233 ),
234 mesh,
236 )
237 );
238
239 volVectorField& Kappa = tKappa.ref();
240 vectorField& KappaCells = Kappa.primitiveFieldRef();
241 const scalarField& TCells = this->T_;
242 const scalarField& pCells = this->p_;
243
244 forAll(KappaCells, celli)
245 {
246 Kappa[celli] =
247 this->cellVolMixture
248 (
249 pCells[celli],
250 TCells[celli],
251 celli
252 ).Kappa(pCells[celli], TCells[celli]);
253 }
254
256
257 forAll(KappaBf, patchi)
258 {
259 vectorField& Kappap = KappaBf[patchi];
260 const scalarField& pT = this->T_.boundaryField()[patchi];
261 const scalarField& pp = this->p_.boundaryField()[patchi];
262
263 forAll(Kappap, facei)
264 {
265 Kappap[facei] =
266 this->patchFaceVolMixture
267 (
268 pp[facei],
269 pT[facei],
270 patchi,
271 facei
272 ).Kappa(pp[facei], pT[facei]);
273 }
274 }
275
276 return tKappa;
277}
278
279
280template<class BasicSolidThermo, class MixtureType>
283(
284 const label patchi
285) const
286{
287 const scalarField& pp = this->p_.boundaryField()[patchi];
288 const scalarField& Tp = this->T_.boundaryField()[patchi];
289 tmp<vectorField> tKappa(new vectorField(pp.size()));
290
291 vectorField& Kappap = tKappa.ref();
292
293 forAll(Tp, facei)
294 {
295 Kappap[facei] =
296 this->patchFaceVolMixture
297 (
298 pp[facei],
299 Tp[facei],
300 patchi,
301 facei
302 ).Kappa(pp[facei], Tp[facei]);
303 }
304
305 return tKappa;
306}
307
308
309// ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Energy for a solid mixture.
Definition: heSolidThermo.H:55
virtual void correct()
Update properties.
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
virtual ~heSolidThermo()
Destructor.
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:56
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: pureMixture.H:66
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
const word dictName("faMeshDefinition")
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333