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 -------------------------------------------------------------------------------
10 License
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 
40 namespace Foam
41 {
42 namespace functionObjects
43 {
44  defineTypeNameAndDebug(energyTransport, 0);
45 
47  (
48  functionObject,
49  energyTransport,
50  dictionary
51  );
52 }
53 }
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
58 Foam::volScalarField& Foam::functionObjects::energyTransport::transportedField()
59 {
60  if (!foundObject<volScalarField>(fieldName_))
61  {
62  auto tfldPtr = tmp<volScalarField>::New
63  (
64  IOobject
65  (
66  fieldName_,
67  mesh_.time().timeName(),
68  mesh_,
72  ),
73  mesh_
74  );
75  store(fieldName_, tfldPtr);
76  }
77 
78  return lookupObjectRef<volScalarField>(fieldName_);
79 }
80 
81 
83 Foam::functionObjects::energyTransport::kappaEff() const
84 {
85  // Incompressible
86  {
87  typedef incompressible::turbulenceModel turbType;
88 
89  const turbType* turbPtr = findObject<turbType>
90  (
92  );
93 
94  if (turbPtr)
95  {
96  return tmp<volScalarField>
97  (
98  new volScalarField
99  (
100  kappa() + Cp()*turbPtr->nut()*rho()/Prt_
101  )
102  );
103  }
104  }
105 
107  << "Turbulence model not found" << exit(FatalError);
108 
109  return nullptr;
110 }
111 
112 
114 Foam::functionObjects::energyTransport::rho() const
115 {
117  (
118  IOobject
119  (
120  "trho",
121  mesh_.time().timeName(),
122  mesh_,
126  ),
127  mesh_,
128  rho_
129  );
130 
131  if (phases_.size())
132  {
133  trho.ref() = lookupObject<volScalarField>(rhoName_);
134  }
135  return trho;
136 }
137 
138 
140 Foam::functionObjects::energyTransport::Cp() const
141 {
142  if (phases_.size())
143  {
144  tmp<volScalarField> tCp(phases_[0]*Cps_[0]);
145 
146  for (label i = 1; i < phases_.size(); i++)
147  {
148  tCp.ref() += phases_[i]*Cps_[i];
149  }
150  return tCp;
151  }
152 
154  (
155  IOobject
156  (
157  "tCp",
158  mesh_.time().timeName(),
159  mesh_,
163  ),
164  mesh_,
165  Cp_
166  );
167 }
168 
169 
171 Foam::functionObjects::energyTransport::kappa() const
172 {
173  if (phases_.size())
174  {
175  tmp<volScalarField> tkappa(phases_[0]*kappas_[0]);
176 
177  for (label i = 1; i < phases_.size(); i++)
178  {
179  tkappa.ref() += phases_[i]*kappas_[i];
180  }
181  return tkappa;
182  }
183 
185  (
186  IOobject
187  (
188  "tkappa",
189  mesh_.time().timeName(),
190  mesh_,
194  ),
195  mesh_,
196  kappa_
197  );
198 }
199 
200 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
201 
202 Foam::functionObjects::energyTransport::energyTransport
203 (
204  const word& name,
205  const Time& runTime,
206  const dictionary& dict
207 )
208 :
210  fieldName_(dict.getOrDefault<word>("field", "T")),
211  phiName_(dict.getOrDefault<word>("phi", "phi")),
212  rhoName_(dict.getOrDefault<word>("rho", "rho")),
213  nCorr_(0),
214  schemesField_("unknown-schemesField"),
215  fvOptions_(mesh_),
216  multiphaseThermo_(dict.subOrEmptyDict("phaseThermos")),
217  Cp_("Cp", dimEnergy/dimMass/dimTemperature, 0, dict),
218  kappa_
219  (
220  "kappa",
222  0,
223  dict
224  ),
225  rho_("rhoInf", dimDensity, 0, dict),
226  Prt_("Prt", dimless, 1, dict),
227  rhoCp_
228  (
229  IOobject
230  (
231  "rhoCp",
232  mesh_.time().timeName(),
233  mesh_,
234  IOobject::NO_READ,
235  IOobject::NO_WRITE,
236  IOobject::NO_REGISTER
237  ),
238  mesh_,
240  )
241 {
242  read(dict);
243 
244  // If the flow is multiphase
245  if (!multiphaseThermo_.empty())
246  {
247  Cps_.setSize(multiphaseThermo_.size());
248  kappas_.setSize(Cps_.size());
249  phaseNames_.setSize(Cps_.size());
250 
251  label phasei = 0;
252  forAllConstIters(multiphaseThermo_, iter)
253  {
254  const word& key = iter().keyword();
255 
256  const dictionary* subDictPtr = multiphaseThermo_.findDict(key);
257 
258  if (!subDictPtr)
259  {
261  << "Found non-dictionary entry " << iter()
262  << " in top-level dictionary " << multiphaseThermo_
263  << exit(FatalError);
264  }
265 
266  const dictionary& dict = *subDictPtr;
267 
268  phaseNames_[phasei] = key;
269 
270  Cps_.set
271  (
272  phasei,
274  (
275  "Cp",
277  dict
278  )
279  );
280 
281  kappas_.set
282  (
283  phasei,
284  new dimensionedScalar //[J/m/s/K]
285  (
286  "kappa",
288  dict
289  )
290  );
291 
292  ++phasei;
293  }
294 
295  phases_.setSize(phaseNames_.size());
296  forAll(phaseNames_, i)
297  {
298  phases_.set
299  (
300  i,
301  mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
302  );
303  }
304 
305  rhoCp_ = rho()*Cp();
306  rhoCp_.oldTime();
307  }
308  else
309  {
310  if (Cp_.value() == 0.0 || kappa_.value() == 0.0)
311  {
313  << " Multiphase thermo dictionary not found and Cp/kappa "
314  << " for single phase are zero. Please entry either"
315  << exit(FatalError);
316  }
317 
318  }
319 
320  // Force creation of transported field so any BCs using it can
321  // look it up
322  volScalarField& s = transportedField();
323  s.correctBoundaryConditions();
324 }
325 
326 
327 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
330 {}
331 
332 
333 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
334 
336 {
338 
339  dict.readIfPresent("phi", phiName_);
340  dict.readIfPresent("rho", rhoName_);
341 
342  schemesField_ = dict.getOrDefault("schemesField", fieldName_);
343 
344  dict.readIfPresent("nCorr", nCorr_);
345 
346  if (dict.found("fvOptions"))
347  {
348  fvOptions_.reset(dict.subDict("fvOptions"));
349  }
350 
351  return true;
352 }
353 
354 
356 {
357  volScalarField& s = transportedField();
358 
359  Log << type() << " execute: " << s.name() << endl;
360 
361  const surfaceScalarField& phi =
362  mesh_.lookupObject<surfaceScalarField>(phiName_);
363 
364  // Calculate the diffusivity
365  const volScalarField kappaEff("kappaEff", this->kappaEff());
366 
367  word divScheme("div(phi," + schemesField_ + ")");
368  word laplacianScheme
369  (
370  "laplacian(kappaEff," + schemesField_ + ")"
371  );
372 
373  // Set under-relaxation coeff
374  scalar relaxCoeff = 0;
375  mesh_.relaxEquation(schemesField_, relaxCoeff);
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(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  (
409  IOobject
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  (
425  fvm::ddt(rhoCp, s)
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(schemesField_);
437  }
438  }
439  else
440  {
442  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
443  << "Dimensions should be " << dimMass/dimTime << " or "
445  }
447  Log << endl;
448 
449  return true;
450 }
451 
452 
454 {
455  return true;
456 }
457 
458 
459 // ************************************************************************* //
Foam::surfaceFields.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
defineTypeNameAndDebug(ObukhovLength, 0)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label phasei
Definition: pEqn.H:27
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
virtual bool read(const dictionary &)
Read the energyTransport data.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< volScalarField > trho
Calculate the matrix for the laplacian of the field.
virtual bool execute()
Calculate the energyTransport.
virtual bool write()
Do nothing.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
kappaEff
Definition: TEqn.H:10
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const surfaceScalarField rhoCpPhi(fvc::interpolate(fluid.Cp()) *rhoPhi)
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
void setSize(const label n)
Alias for resize()
Definition: List.H:316
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
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...
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
Calculate the matrix for the first temporal derivative.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
rhoCp
Definition: TEqn.H:3
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
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
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
const volScalarField & Cp
Definition: EEqn.H:7
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1101
const dimensionSet dimEnergy
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.
const dimensionSet dimDensity
Calculate the matrix for the divergence of the given field and flux.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
Nothing to be read.
#define Log
Definition: PDRblock.C:28
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual bool read(const dictionary &dict)
Read optional controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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))
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
Do not request registration (bool: false)
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127