Loading...
Searching...
No Matches
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-2023 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"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35namespace Foam
36{
37namespace functionObjects
38{
40
42 (
46 );
47}
48}
49
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
53Foam::volScalarField& Foam::functionObjects::energyTransport::transportedField()
54{
55 if (!foundObject<volScalarField>(fieldName_))
56 {
57 auto tfldPtr = tmp<volScalarField>::New
58 (
60 (
61 fieldName_,
63 mesh_,
67 ),
68 mesh_
69 );
70 store(fieldName_, tfldPtr);
71 }
72
73 return lookupObjectRef<volScalarField>(fieldName_);
74}
75
76
77Foam::tmp<Foam::volScalarField>
78Foam::functionObjects::energyTransport::kappaEff() const
79{
80 // Incompressible
81 {
82 typedef incompressible::turbulenceModel turbType;
83
84 const turbType* turbPtr = findObject<turbType>
85 (
87 );
88
89 if (turbPtr)
90 {
92 (
93 "kappaEff",
95 (
96 kappa() + Cp()*turbPtr->nut()*rho()/Prt_
97 )
98 );
99 }
100 }
101
103 << "Turbulence model not found" << exit(FatalError);
104
105 return nullptr;
106}
107
108
109Foam::tmp<Foam::volScalarField>
110Foam::functionObjects::energyTransport::rho() const
111{
113 (
114 "trho",
116 mesh_,
117 rho_
118 );
119
120 if (phases_.size())
121 {
122 trho.ref() = lookupObject<volScalarField>(rhoName_);
123 }
124 return trho;
125}
126
127
128Foam::tmp<Foam::volScalarField>
129Foam::functionObjects::energyTransport::Cp() const
130{
131 if (phases_.size())
132 {
133 tmp<volScalarField> tCp(phases_[0]*Cps_[0]);
134
135 for (label i = 1; i < phases_.size(); i++)
136 {
137 tCp.ref() += phases_[i]*Cps_[i];
138 }
139 return tCp;
140 }
141
143 (
144 "tCp",
146 mesh_,
147 Cp_
148 );
149}
150
151
152Foam::tmp<Foam::volScalarField>
153Foam::functionObjects::energyTransport::kappa() const
154{
155 if (phases_.size())
156 {
157 tmp<volScalarField> tkappa(phases_[0]*kappas_[0]);
158
159 for (label i = 1; i < phases_.size(); i++)
160 {
161 tkappa.ref() += phases_[i]*kappas_[i];
162 }
163 return tkappa;
164 }
165
167 (
168 "tkappa",
170 mesh_,
171 kappa_
172 );
173}
174
175
176// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
177
179(
180 const word& name,
181 const Time& runTime,
182 const dictionary& dict
183)
184:
186 rhoCp_
187 (
189 (
190 "rhoCp",
191 mesh_.time().timeName(),
192 mesh_,
193 IOobject::NO_READ,
194 IOobject::NO_WRITE,
195 IOobject::NO_REGISTER
196 ),
197 mesh_,
199 ),
200 fvOptions_(mesh_),
201 multiphaseThermo_(dict.subOrEmptyDict("phaseThermos")),
202 Cp_("Cp", dimEnergy/dimMass/dimTemperature, 0, dict),
203 kappa_
204 (
205 "kappa",
207 0,
208 dict
209 ),
210 rho_("rhoInf", dimDensity, 0, dict),
211 Prt_("Prt", dimless, 1, dict),
212 fieldName_(dict.getOrDefault<word>("field", "T")),
213 schemesField_("unknown-schemesField"),
214 phiName_(dict.getOrDefault<word>("phi", "phi")),
215 rhoName_(dict.getOrDefault<word>("rho", "rho")),
216 tol_(1),
217 nCorr_(0)
218{
219 read(dict);
220
221 // If the flow is multiphase
222 if (!multiphaseThermo_.empty())
223 {
224 Cps_.setSize(multiphaseThermo_.size());
225 kappas_.setSize(Cps_.size());
226 phaseNames_.setSize(Cps_.size());
227
228 label phasei = 0;
229 forAllConstIters(multiphaseThermo_, iter)
230 {
231 const word& key = iter().keyword();
232
233 const dictionary* subDictPtr = multiphaseThermo_.findDict(key);
234
235 if (!subDictPtr)
236 {
238 << "Found non-dictionary entry " << iter()
239 << " in top-level dictionary " << multiphaseThermo_
240 << exit(FatalError);
241 }
242
243 const dictionary& dict = *subDictPtr;
244
245 phaseNames_[phasei] = key;
246
247 Cps_.set
248 (
249 phasei,
251 (
252 "Cp",
254 dict
255 )
256 );
257
258 kappas_.set
259 (
260 phasei,
261 new dimensionedScalar //[J/m/s/K]
262 (
263 "kappa",
265 dict
266 )
267 );
268
269 ++phasei;
270 }
271
272 phases_.setSize(phaseNames_.size());
273 forAll(phaseNames_, i)
274 {
275 phases_.set
276 (
277 i,
278 mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
279 );
280 }
281
282 rhoCp_ = rho()*Cp();
283 rhoCp_.oldTime();
284 }
285 else
286 {
287 if (Cp_.value() == 0.0 || kappa_.value() == 0.0)
288 {
290 << " Multiphase thermo dictionary not found and Cp/kappa "
291 << " for single phase are zero. Please entry either"
292 << exit(FatalError);
293 }
294
295 }
296
297 // Force creation of transported field so any BCs using it can
298 // look it up
299 volScalarField& s = transportedField();
300 s.correctBoundaryConditions();
301}
302
303
304// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
305
307{
309 {
310 return false;
311 }
312
313 dict.readIfPresent("phi", phiName_);
314 dict.readIfPresent("rho", rhoName_);
315
316 schemesField_ = dict.getOrDefault("schemesField", fieldName_);
317
318 dict.readIfPresent("nCorr", nCorr_);
319 dict.readIfPresent("tolerance", tol_);
320
321 if (dict.found("fvOptions"))
322 {
323 fvOptions_.reset(dict.subDict("fvOptions"));
324 }
325
326 return true;
327}
328
329
331{
332 volScalarField& s = transportedField();
333
334 Log << type() << " execute: " << s.name() << endl;
335
336 const surfaceScalarField& phi =
337 mesh_.lookupObject<surfaceScalarField>(phiName_);
338
339 // Calculate the diffusivity
340 const volScalarField kappaEff("kappaEff", this->kappaEff());
341
342 word divScheme("div(phi," + schemesField_ + ")");
343 word laplacianScheme
344 (
345 "laplacian(kappaEff," + schemesField_ + ")"
346 );
347
348 // Set under-relaxation coeff
349 scalar relaxCoeff = 0;
350 mesh_.relaxEquation(schemesField_, relaxCoeff);
351
352 // Convergence monitor parameters
353 bool converged = false;
354 int iter = 0;
355
356 if (phi.dimensions() == dimMass/dimTime)
357 {
358 rhoCp_ = rho()*Cp();
360
361 for (int i = 0; i <= nCorr_; ++i)
362 {
363 fvScalarMatrix sEqn
364 (
365 fvm::ddt(rhoCp_, s)
366 + fvm::div(rhoCpPhi, s, divScheme)
367 - fvm::Sp(fvc::ddt(rhoCp_) + fvc::div(rhoCpPhi), s)
368 - fvm::laplacian(kappaEff, s, laplacianScheme)
369 ==
370 fvOptions_(rhoCp_, s)
371 );
372
373 sEqn.relax(relaxCoeff);
374
375 fvOptions_.constrain(sEqn);
376
377 ++iter;
378 converged = (sEqn.solve(schemesField_).initialResidual() < tol_);
379 if (converged) break;
380 }
381 }
382 else if (phi.dimensions() == dimVolume/dimTime)
383 {
384 dimensionedScalar rhoCp(rho_*Cp_);
385
386 const surfaceScalarField CpPhi(rhoCp*phi);
387
388 auto trhoCp = volScalarField::New
389 (
390 "trhoCp",
392 mesh_,
393 rhoCp
394 );
395
396 for (int i = 0; i <= nCorr_; ++i)
397 {
398 fvScalarMatrix sEqn
399 (
401 + fvm::div(CpPhi, s, divScheme)
402 - fvm::laplacian(kappaEff, s, laplacianScheme)
403 ==
404 fvOptions_(trhoCp.ref(), s)
405 );
406
407 sEqn.relax(relaxCoeff);
408
409 fvOptions_.constrain(sEqn);
410
411 ++iter;
412 converged = (sEqn.solve(schemesField_).initialResidual() < tol_);
413 if (converged) break;
414 }
415 }
416 else
417 {
419 << "Incompatible dimensions for phi: " << phi.dimensions() << nl
420 << "Dimensions should be " << dimMass/dimTime << " or "
422 }
423
424 if (converged)
425 {
426 Log << type() << ": " << name() << ": "
427 << s.name() << " is converged." << nl
428 << tab << "initial-residual tolerance: " << tol_ << nl
429 << tab << "outer iteration: " << iter << nl;
430 }
432 Log << endl;
433
434 return true;
435}
436
437
439{
440 Log << type() << " write: " << name() << nl
441 << tab << fieldName_ << nl
442 << endl;
443
444 volScalarField& s = transportedField();
445
446 s.write();
447
448 return true;
449}
450
451
452// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
bool empty() const noexcept
True if the list is empty.
Definition DLListBase.H:189
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void setSize(label n)
Alias for resize().
Definition List.H:536
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Computes the simplified energy transport equation in single-phase or two-phase flow,...
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function object output.
energyTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
virtual bool read(const dictionary &)
Read the function-object dictionary.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
bool foundObject(const word &fieldName) const
Find object (eg, a field) in the (sub) objectRegistry.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
ObjectType & lookupObjectRef(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
const Time & time() const
Return time database.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
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:600
auto & name
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].
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
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:42
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
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)
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
tmp< volScalarField > trho
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235