energyTransport.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) 2017-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "energyTransport.H"
29#include "surfaceFields.H"
30#include "fvmDdt.H"
31#include "fvmDiv.H"
32#include "fvmLaplacian.H"
33#include "fvmSup.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42namespace functionObjects
43{
45
47 (
51 );
52}
53}
54
55
56// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57
58Foam::volScalarField& Foam::functionObjects::energyTransport::transportedField()
59{
60 if (!foundObject<volScalarField>(fieldName_))
61 {
62 auto tfldPtr = tmp<volScalarField>::New
63 (
64 IOobject
65 (
66 fieldName_,
68 mesh_,
71 ),
72 mesh_
73 );
74 store(fieldName_, tfldPtr);
75 }
76
77 return lookupObjectRef<volScalarField>(fieldName_);
78}
79
80
82Foam::functionObjects::energyTransport::kappaEff() const
83{
84 // Incompressible
85 {
86 typedef incompressible::turbulenceModel turbType;
87
88 const turbType* turbPtr = findObject<turbType>
89 (
91 );
92
93 if (turbPtr)
94 {
95 return tmp<volScalarField>
96 (
98 (
99 kappa() + Cp()*turbPtr->nut()*rho()/Prt_
100 )
101 );
102 }
103 }
104
106 << "Turbulence model not found" << exit(FatalError);
107
108 return nullptr;
109}
110
111
113Foam::functionObjects::energyTransport::rho() const
114{
116 (
117 IOobject
118 (
119 "trho",
120 mesh_.time().timeName(),
121 mesh_,
124 false
125 ),
126 mesh_,
127 rho_
128 );
129
130 if (phases_.size())
131 {
132 trho.ref() = lookupObject<volScalarField>(rhoName_);
133 }
134 return trho;
135}
136
137
139Foam::functionObjects::energyTransport::Cp() const
140{
141 if (phases_.size())
142 {
143 tmp<volScalarField> tCp(phases_[0]*Cps_[0]);
144
145 for (label i = 1; i < phases_.size(); i++)
146 {
147 tCp.ref() += phases_[i]*Cps_[i];
148 }
149 return tCp;
150 }
151
153 (
154 IOobject
155 (
156 "tCp",
157 mesh_.time().timeName(),
158 mesh_,
161 false
162 ),
163 mesh_,
164 Cp_
165 );
166}
167
168
170Foam::functionObjects::energyTransport::kappa() const
171{
172 if (phases_.size())
173 {
174 tmp<volScalarField> tkappa(phases_[0]*kappas_[0]);
175
176 for (label i = 1; i < phases_.size(); i++)
177 {
178 tkappa.ref() += phases_[i]*kappas_[i];
179 }
180 return tkappa;
181 }
182
184 (
185 IOobject
186 (
187 "tkappa",
188 mesh_.time().timeName(),
189 mesh_,
192 false
193 ),
194 mesh_,
195 kappa_
196 );
197}
198
199// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
200
202(
203 const word& name,
204 const Time& runTime,
205 const dictionary& dict
206)
207:
209 fieldName_(dict.getOrDefault<word>("field", "T")),
210 phiName_(dict.getOrDefault<word>("phi", "phi")),
211 rhoName_(dict.getOrDefault<word>("rho", "rho")),
212 nCorr_(0),
213 schemesField_("unknown-schemesField"),
214 fvOptions_(mesh_),
215 multiphaseThermo_(dict.subOrEmptyDict("phaseThermos")),
216 Cp_("Cp", dimEnergy/dimMass/dimTemperature, 0, dict),
217 kappa_
218 (
219 "kappa",
221 0,
222 dict
223 ),
224 rho_("rhoInf", dimDensity, 0, dict),
225 Prt_("Prt", dimless, 1, dict),
226 rhoCp_
227 (
229 (
230 "rhoCp",
231 mesh_.time().timeName(),
232 mesh_,
233 IOobject::NO_READ,
234 IOobject::NO_WRITE,
235 false
236 ),
237 mesh_,
239 )
240{
241 read(dict);
242
243 // If the flow is multiphase
244 if (!multiphaseThermo_.empty())
245 {
246 Cps_.setSize(multiphaseThermo_.size());
247 kappas_.setSize(Cps_.size());
248 phaseNames_.setSize(Cps_.size());
249
250 label phasei = 0;
251 forAllConstIters(multiphaseThermo_, iter)
252 {
253 const word& key = iter().keyword();
254
255 if (!multiphaseThermo_.isDict(key))
256 {
258 << "Found non-dictionary entry " << iter()
259 << " in top-level dictionary " << multiphaseThermo_
260 << exit(FatalError);
261 }
262
263 const dictionary& dict = multiphaseThermo_.subDict(key);
264
265 phaseNames_[phasei] = key;
266
267 Cps_.set
268 (
269 phasei,
271 (
272 "Cp",
274 dict
275 )
276 );
277
278 kappas_.set
279 (
280 phasei,
281 new dimensionedScalar //[J/m/s/K]
282 (
283 "kappa",
285 dict
286 )
287 );
288
289 ++phasei;
290 }
291
292 phases_.setSize(phaseNames_.size());
293 forAll(phaseNames_, i)
294 {
295 phases_.set
296 (
297 i,
298 mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
299 );
300 }
301
302 rhoCp_ = rho()*Cp();
303 rhoCp_.oldTime();
304 }
305 else
306 {
307 if (Cp_.value() == 0.0 || kappa_.value() == 0.0)
308 {
310 << " Multiphase thermo dictionary not found and Cp/kappa "
311 << " for single phase are zero. Please entry either"
312 << exit(FatalError);
313 }
314
315 }
316
317 // Force creation of transported field so any BCs using it can
318 // look it up
319 volScalarField& s = transportedField();
320 s.correctBoundaryConditions();
321}
322
323
324// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
325
327{}
328
329
330// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
331
333{
335
336 dict.readIfPresent("phi", phiName_);
337 dict.readIfPresent("rho", rhoName_);
338
339 schemesField_ = dict.getOrDefault("schemesField", fieldName_);
340
341 dict.readIfPresent("nCorr", nCorr_);
342
343 if (dict.found("fvOptions"))
344 {
345 fvOptions_.reset(dict.subDict("fvOptions"));
346 }
347
348 return true;
349}
350
351
353{
354 volScalarField& s = transportedField();
355
356 Log << type() << " execute: " << s.name() << endl;
357
358 const surfaceScalarField& phi =
359 mesh_.lookupObject<surfaceScalarField>(phiName_);
360
361 // Calculate the diffusivity
362 const volScalarField kappaEff("kappaEff", this->kappaEff());
363
364 word divScheme("div(phi," + schemesField_ + ")");
365 word laplacianScheme
366 (
367 "laplacian(kappaEff," + schemesField_ + ")"
368 );
369
370 // Set under-relaxation coeff
371 scalar relaxCoeff = 0.0;
372 if (mesh_.relaxEquation(schemesField_))
373 {
374 relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
375 }
376
377 if (phi.dimensions() == dimMass/dimTime)
378 {
379 rhoCp_ = rho()*Cp();
381
382 for (label i = 0; i <= nCorr_; i++)
383 {
384 fvScalarMatrix sEqn
385 (
386 fvm::ddt(rhoCp_, s)
387 + fvm::div(rhoCpPhi, s, divScheme)
388 - fvm::Sp(fvc::ddt(rhoCp_) + fvc::div(rhoCpPhi), s)
389 - fvm::laplacian(kappaEff, s, laplacianScheme)
390 ==
391 fvOptions_(rhoCp_, s)
392 );
393
394 sEqn.relax(relaxCoeff);
395
396 fvOptions_.constrain(sEqn);
397
398 sEqn.solve(mesh_.solverDict(schemesField_));
399 }
400 }
401 else if (phi.dimensions() == dimVolume/dimTime)
402 {
403 dimensionedScalar rhoCp(rho_*Cp_);
404
405 const surfaceScalarField CpPhi(rhoCp*phi);
406
407 auto trhoCp = tmp<volScalarField>::New
408 (
410 (
411 "trhoCp",
412 mesh_.time().timeName(),
413 mesh_,
416 ),
417 mesh_,
418 rhoCp
419 );
420
421 for (label i = 0; i <= nCorr_; i++)
422 {
423 fvScalarMatrix sEqn
424 (
426 + fvm::div(CpPhi, s, divScheme)
427 - fvm::laplacian(kappaEff, s, laplacianScheme)
428 ==
429 fvOptions_(trhoCp.ref(), s)
430 );
431
432 sEqn.relax(relaxCoeff);
433
434 fvOptions_.constrain(sEqn);
435
436 sEqn.solve(mesh_.solverDict(schemesField_));
437 }
438 }
439 else
440 {
442 << "Incompatible dimensions for phi: " << phi.dimensions() << nl
443 << "Dimensions should be " << dimMass/dimTime << " or "
445 }
446
447 Log << endl;
448
449 return true;
450}
451
452
454{
455 return true;
456}
457
458
459// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void setSize(const label n)
Alias for resize()
Definition: List.H:218
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
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
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
bool isDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Check if entry is found and is a sub-dictionary.
Definition: dictionaryI.H:147
const Type & value() const
Return const reference to value.
Abstract base-class for Time/database function objects.
Evolves a simplified energy transport equation for incompressible flows. It takes into account the in...
virtual bool execute()
Calculate the energyTransport.
virtual bool read(const dictionary &)
Read the energyTransport data.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Type * getObjectPtr(const word &name, const bool recursive=false) const
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
const volScalarField & Cp
Definition: EEqn.H:7
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
word timeName
Definition: getTimeIndex.H:3
label phasei
Definition: pEqn.H:27
const surfaceScalarField rhoCpPhi(fvc::interpolate(fluid.Cp()) *rhoPhi)
rhoCp
Definition: TEqn.H:3
kappaEff
Definition: TEqn.H:10
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
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
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
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
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
tmp< volScalarField > trho
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Foam::surfaceFields.