BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / thermophysicalModels / basic / basicThermo / basicThermo.C
blob7b0d81aea78eb41e44baba35b885f283a451212f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "basicThermo.H"
27 #include "fvMesh.H"
28 #include "HashTable.H"
29 #include "zeroGradientFvPatchFields.H"
30 #include "fixedEnthalpyFvPatchScalarField.H"
31 #include "gradientEnthalpyFvPatchScalarField.H"
32 #include "mixedEnthalpyFvPatchScalarField.H"
33 #include "fixedInternalEnergyFvPatchScalarField.H"
34 #include "gradientInternalEnergyFvPatchScalarField.H"
35 #include "mixedInternalEnergyFvPatchScalarField.H"
37 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
39 namespace Foam
41     defineTypeNameAndDebug(basicThermo, 0);
44 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
46 Foam::wordList Foam::basicThermo::hBoundaryTypes()
48     const volScalarField::GeometricBoundaryField& tbf = T_.boundaryField();
50     wordList hbt = tbf.types();
52     forAll(tbf, patchi)
53     {
54         if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
55         {
56             hbt[patchi] = fixedEnthalpyFvPatchScalarField::typeName;
57         }
58         else if
59         (
60             isA<zeroGradientFvPatchScalarField>(tbf[patchi])
61          || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
62         )
63         {
64             hbt[patchi] = gradientEnthalpyFvPatchScalarField::typeName;
65         }
66         else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
67         {
68             hbt[patchi] = mixedEnthalpyFvPatchScalarField::typeName;
69         }
70     }
72     return hbt;
76 void Foam::basicThermo::hBoundaryCorrection(volScalarField& h)
78     volScalarField::GeometricBoundaryField& hbf = h.boundaryField();
80     forAll(hbf, patchi)
81     {
82         if (isA<gradientEnthalpyFvPatchScalarField>(hbf[patchi]))
83         {
84             refCast<gradientEnthalpyFvPatchScalarField>(hbf[patchi]).gradient()
85                 = hbf[patchi].fvPatchField::snGrad();
86         }
87         else if (isA<mixedEnthalpyFvPatchScalarField>(hbf[patchi]))
88         {
89             refCast<mixedEnthalpyFvPatchScalarField>(hbf[patchi]).refGrad()
90                 = hbf[patchi].fvPatchField::snGrad();
91         }
92     }
96 Foam::wordList Foam::basicThermo::eBoundaryTypes()
98     const volScalarField::GeometricBoundaryField& tbf = T_.boundaryField();
100     wordList ebt = tbf.types();
102     forAll(tbf, patchi)
103     {
104         if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
105         {
106             ebt[patchi] = fixedInternalEnergyFvPatchScalarField::typeName;
107         }
108         else if
109         (
110             isA<zeroGradientFvPatchScalarField>(tbf[patchi])
111          || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
112         )
113         {
114             ebt[patchi] = gradientInternalEnergyFvPatchScalarField::typeName;
115         }
116         else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
117         {
118             ebt[patchi] = mixedInternalEnergyFvPatchScalarField::typeName;
119         }
120     }
122     return ebt;
126 void Foam::basicThermo::eBoundaryCorrection(volScalarField& e)
128     volScalarField::GeometricBoundaryField& ebf = e.boundaryField();
130     forAll(ebf, patchi)
131     {
132         if (isA<gradientInternalEnergyFvPatchScalarField>(ebf[patchi]))
133         {
134             refCast<gradientInternalEnergyFvPatchScalarField>(ebf[patchi])
135                 .gradient() = ebf[patchi].fvPatchField::snGrad();
136         }
137         else if (isA<mixedInternalEnergyFvPatchScalarField>(ebf[patchi]))
138         {
139             refCast<mixedInternalEnergyFvPatchScalarField>(ebf[patchi])
140                 .refGrad() = ebf[patchi].fvPatchField::snGrad();
141         }
142     }
145 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
147 Foam::basicThermo::basicThermo(const fvMesh& mesh)
149     IOdictionary
150     (
151         IOobject
152         (
153             "thermophysicalProperties",
154             mesh.time().constant(),
155             mesh,
156             IOobject::MUST_READ_IF_MODIFIED,
157             IOobject::NO_WRITE
158         )
159     ),
161     p_
162     (
163         IOobject
164         (
165             "p",
166             mesh.time().timeName(),
167             mesh,
168             IOobject::MUST_READ,
169             IOobject::AUTO_WRITE
170         ),
171         mesh
172     ),
174     psi_
175     (
176         IOobject
177         (
178             "psi",
179             mesh.time().timeName(),
180             mesh,
181             IOobject::NO_READ,
182             IOobject::NO_WRITE
183         ),
184         mesh,
185         dimensionSet(0, -2, 2, 0, 0)
186     ),
188     T_
189     (
190         IOobject
191         (
192             "T",
193             mesh.time().timeName(),
194             mesh,
195             IOobject::MUST_READ,
196             IOobject::AUTO_WRITE
197         ),
198         mesh
199     ),
201     mu_
202     (
203         IOobject
204         (
205             "mu",
206             mesh.time().timeName(),
207             mesh,
208             IOobject::NO_READ,
209             IOobject::NO_WRITE
210         ),
211         mesh,
212         dimensionSet(1, -1, -1, 0, 0)
213     ),
215     alpha_
216     (
217         IOobject
218         (
219             "alpha",
220             mesh.time().timeName(),
221             mesh,
222             IOobject::NO_READ,
223             IOobject::NO_WRITE
224         ),
225         mesh,
226         dimensionSet(1, -1, -1, 0, 0)
227     )
231 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
233 Foam::basicThermo::~basicThermo()
237 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
239 Foam::volScalarField& Foam::basicThermo::p()
241     return p_;
245 const Foam::volScalarField& Foam::basicThermo::p() const
247     return p_;
251 const Foam::volScalarField& Foam::basicThermo::psi() const
253     return psi_;
257 Foam::volScalarField& Foam::basicThermo::h()
259     notImplemented("basicThermo::h()");
260     return const_cast<volScalarField&>(volScalarField::null());
264 const Foam::volScalarField& Foam::basicThermo::h() const
266     notImplemented("basicThermo::h() const");
267     return volScalarField::null();
271 Foam::tmp<Foam::scalarField> Foam::basicThermo::h
273     const scalarField& T,
274     const labelList& cells
275 ) const
277     notImplemented
278     (
279         "basicThermo::h"
280         "(const scalarField& T, const labelList& cells) const"
281     );
282     return tmp<scalarField>(NULL);
286 Foam::tmp<Foam::scalarField> Foam::basicThermo::h
288     const scalarField& T,
289     const label patchi
290 ) const
292     notImplemented
293     (
294         "basicThermo::h"
295         "(const scalarField& T, const label patchi) const"
296     );
297     return tmp<scalarField>(NULL);
301 Foam::volScalarField& Foam::basicThermo::hs()
303     notImplemented("basicThermo::hs()");
304     return const_cast<volScalarField&>(volScalarField::null());
308 const Foam::volScalarField& Foam::basicThermo::hs() const
310     notImplemented("basicThermo::hs() const");
311     return volScalarField::null();
315 Foam::tmp<Foam::scalarField> Foam::basicThermo::hs
317     const scalarField& T,
318     const labelList& cells
319 ) const
321     notImplemented
322     (
323         "basicThermo::hs"
324         "(const scalarField& T, const labelList& cells) const"
325     );
326     return tmp<scalarField>(NULL);
330 Foam::tmp<Foam::scalarField> Foam::basicThermo::hs
332     const scalarField& T,
333     const label patchi
334 ) const
336     notImplemented
337     (
338         "basicThermo::hs"
339         "(const scalarField& T, const label patchi) const"
340     );
341     return tmp<scalarField>(NULL);
345 Foam::tmp<Foam::volScalarField> Foam::basicThermo::hc() const
347     notImplemented("basicThermo::hc()");
348     return volScalarField::null();
352 Foam::volScalarField& Foam::basicThermo::e()
354     notImplemented("basicThermo::e()");
355     return const_cast<volScalarField&>(volScalarField::null());
359 const Foam::volScalarField& Foam::basicThermo::e() const
361     notImplemented("basicThermo::e()");
362     return volScalarField::null();
366 Foam::tmp<Foam::scalarField> Foam::basicThermo::e
368     const scalarField& T,
369     const labelList& cells
370 ) const
372     notImplemented
373     (
374         "basicThermo::e"
375         "(const scalarField& T, const labelList& cells) const"
376     );
377     return tmp<scalarField>(NULL);
381 Foam::tmp<Foam::scalarField> Foam::basicThermo::e
383     const scalarField& T,
384     const label patchi
385 ) const
387     notImplemented
388     (
389         "basicThermo::e"
390         "(const scalarField& T, const label patchi) const"
391     );
392     return tmp<scalarField>(NULL);
396 const Foam::volScalarField& Foam::basicThermo::T() const
398     return T_;
402 Foam::tmp<Foam::scalarField> Foam::basicThermo::Cp
404     const scalarField& T,
405     const label patchi
406 ) const
408     notImplemented
409     (
410         "basicThermo::Cp"
411         "(const scalarField& T, const label patchi) const"
412     );
413     return tmp<scalarField>(NULL);
417 Foam::tmp<Foam::volScalarField> Foam::basicThermo::Cp() const
419     notImplemented("basicThermo::Cp() const");
420     return volScalarField::null();
424 Foam::tmp<Foam::scalarField> Foam::basicThermo::Cv
426     const scalarField& T,
427     const label patchi
428 ) const
430     notImplemented
431     (
432         "basicThermo::Cv"
433         "(const scalarField& T, const label patchi) const"
434     );
435     return tmp<scalarField>(NULL);
439 Foam::tmp<Foam::volScalarField> Foam::basicThermo::Cv() const
441     notImplemented("basicThermo::Cv() const");
442     return volScalarField::null();
446 const Foam::volScalarField& Foam::basicThermo::mu() const
448     return mu_;
452 const Foam::volScalarField& Foam::basicThermo::alpha() const
454     return alpha_;
458 bool Foam::basicThermo::read()
460     return regIOobject::read();
464 // ************************************************************************* //