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-2022 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_,
71  ),
72  mesh_
73  );
74  store(fieldName_, tfldPtr);
75  }
76 
77  return lookupObjectRef<volScalarField>(fieldName_);
78 }
79 
80 
82 Foam::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  (
97  new volScalarField
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 
113 Foam::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 
139 Foam::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 
170 Foam::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 
201 Foam::functionObjects::energyTransport::energyTransport
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  (
228  IOobject
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  const dictionary* subDictPtr = multiphaseThermo_.findDict(key);
256 
257  if (!subDictPtr)
258  {
260  << "Found non-dictionary entry " << iter()
261  << " in top-level dictionary " << multiphaseThermo_
262  << exit(FatalError);
263  }
264 
265  const dictionary& dict = *subDictPtr;
266 
267  phaseNames_[phasei] = key;
268 
269  Cps_.set
270  (
271  phasei,
273  (
274  "Cp",
276  dict
277  )
278  );
279 
280  kappas_.set
281  (
282  phasei,
283  new dimensionedScalar //[J/m/s/K]
284  (
285  "kappa",
287  dict
288  )
289  );
290 
291  ++phasei;
292  }
293 
294  phases_.setSize(phaseNames_.size());
295  forAll(phaseNames_, i)
296  {
297  phases_.set
298  (
299  i,
300  mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
301  );
302  }
303 
304  rhoCp_ = rho()*Cp();
305  rhoCp_.oldTime();
306  }
307  else
308  {
309  if (Cp_.value() == 0.0 || kappa_.value() == 0.0)
310  {
312  << " Multiphase thermo dictionary not found and Cp/kappa "
313  << " for single phase are zero. Please entry either"
314  << exit(FatalError);
315  }
316 
317  }
318 
319  // Force creation of transported field so any BCs using it can
320  // look it up
321  volScalarField& s = transportedField();
322  s.correctBoundaryConditions();
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
329 {}
330 
331 
332 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
333 
335 {
337 
338  dict.readIfPresent("phi", phiName_);
339  dict.readIfPresent("rho", rhoName_);
340 
341  schemesField_ = dict.getOrDefault("schemesField", fieldName_);
342 
343  dict.readIfPresent("nCorr", nCorr_);
344 
345  if (dict.found("fvOptions"))
346  {
347  fvOptions_.reset(dict.subDict("fvOptions"));
348  }
349 
350  return true;
351 }
352 
353 
355 {
356  volScalarField& s = transportedField();
357 
358  Log << type() << " execute: " << s.name() << endl;
359 
360  const surfaceScalarField& phi =
361  mesh_.lookupObject<surfaceScalarField>(phiName_);
362 
363  // Calculate the diffusivity
364  const volScalarField kappaEff("kappaEff", this->kappaEff());
365 
366  word divScheme("div(phi," + schemesField_ + ")");
367  word laplacianScheme
368  (
369  "laplacian(kappaEff," + schemesField_ + ")"
370  );
371 
372  // Set under-relaxation coeff
373  scalar relaxCoeff = 0.0;
374  if (mesh_.relaxEquation(schemesField_))
375  {
376  relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
377  }
378 
379  if (phi.dimensions() == dimMass/dimTime)
380  {
381  rhoCp_ = rho()*Cp();
383 
384  for (label i = 0; i <= nCorr_; i++)
385  {
386  fvScalarMatrix sEqn
387  (
388  fvm::ddt(rhoCp_, s)
389  + fvm::div(rhoCpPhi, s, divScheme)
390  - fvm::Sp(fvc::ddt(rhoCp_) + fvc::div(rhoCpPhi), s)
391  - fvm::laplacian(kappaEff, s, laplacianScheme)
392  ==
393  fvOptions_(rhoCp_, s)
394  );
395 
396  sEqn.relax(relaxCoeff);
397 
398  fvOptions_.constrain(sEqn);
399 
400  sEqn.solve(mesh_.solverDict(schemesField_));
401  }
402  }
403  else if (phi.dimensions() == dimVolume/dimTime)
404  {
405  dimensionedScalar rhoCp(rho_*Cp_);
406 
407  const surfaceScalarField CpPhi(rhoCp*phi);
408 
409  auto trhoCp = tmp<volScalarField>::New
410  (
411  IOobject
412  (
413  "trhoCp",
414  mesh_.time().timeName(),
415  mesh_,
418  ),
419  mesh_,
420  rhoCp
421  );
422 
423  for (label i = 0; i <= nCorr_; i++)
424  {
425  fvScalarMatrix sEqn
426  (
427  fvm::ddt(rhoCp, s)
428  + fvm::div(CpPhi, s, divScheme)
429  - fvm::laplacian(kappaEff, s, laplacianScheme)
430  ==
431  fvOptions_(trhoCp.ref(), s)
432  );
433 
434  sEqn.relax(relaxCoeff);
435 
436  fvOptions_.constrain(sEqn);
437 
438  sEqn.solve(mesh_.solverDict(schemesField_));
439  }
440  }
441  else
442  {
444  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
445  << "Dimensions should be " << dimMass/dimTime << " or "
447  }
449  Log << endl;
450 
451  return true;
452 }
453 
454 
456 {
457  return true;
458 }
459 
460 
461 // ************************************************************************* //
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:118
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:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
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:49
engineTime & runTime
virtual bool read(const dictionary &)
Read the energyTransport data.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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:361
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:413
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:84
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:752
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
void setSize(const label n)
Alias for resize()
Definition: List.H:289
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
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:203
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
rhoCp
Definition: TEqn.H:3
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
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)
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:166
const fvMesh & mesh_
Reference to the fvMesh.
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:120
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157