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 volScalarField::New
97  (
98  "kappaEff",
100  (
101  kappa() + Cp()*turbPtr->nut()*rho()/Prt_
102  )
103  );
104  }
105  }
106 
108  << "Turbulence model not found" << exit(FatalError);
109 
110  return nullptr;
111 }
112 
113 
115 Foam::functionObjects::energyTransport::rho() const
116 {
118  (
119  "trho",
121  mesh_,
122  rho_
123  );
124 
125  if (phases_.size())
126  {
127  trho.ref() = lookupObject<volScalarField>(rhoName_);
128  }
129  return trho;
130 }
131 
132 
134 Foam::functionObjects::energyTransport::Cp() const
135 {
136  if (phases_.size())
137  {
138  tmp<volScalarField> tCp(phases_[0]*Cps_[0]);
139 
140  for (label i = 1; i < phases_.size(); i++)
141  {
142  tCp.ref() += phases_[i]*Cps_[i];
143  }
144  return tCp;
145  }
146 
147  return volScalarField::New
148  (
149  "tCp",
151  mesh_,
152  Cp_
153  );
154 }
155 
156 
158 Foam::functionObjects::energyTransport::kappa() const
159 {
160  if (phases_.size())
161  {
162  tmp<volScalarField> tkappa(phases_[0]*kappas_[0]);
163 
164  for (label i = 1; i < phases_.size(); i++)
165  {
166  tkappa.ref() += phases_[i]*kappas_[i];
167  }
168  return tkappa;
169  }
170 
171  return volScalarField::New
172  (
173  "tkappa",
175  mesh_,
176  kappa_
177  );
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
182 
183 Foam::functionObjects::energyTransport::energyTransport
184 (
185  const word& name,
186  const Time& runTime,
187  const dictionary& dict
188 )
189 :
191  fieldName_(dict.getOrDefault<word>("field", "T")),
192  phiName_(dict.getOrDefault<word>("phi", "phi")),
193  rhoName_(dict.getOrDefault<word>("rho", "rho")),
194  nCorr_(0),
195  schemesField_("unknown-schemesField"),
196  fvOptions_(mesh_),
197  multiphaseThermo_(dict.subOrEmptyDict("phaseThermos")),
198  Cp_("Cp", dimEnergy/dimMass/dimTemperature, 0, dict),
199  kappa_
200  (
201  "kappa",
203  0,
204  dict
205  ),
206  rho_("rhoInf", dimDensity, 0, dict),
207  Prt_("Prt", dimless, 1, dict),
208  rhoCp_
209  (
210  IOobject
211  (
212  "rhoCp",
213  mesh_.time().timeName(),
214  mesh_,
215  IOobject::NO_READ,
216  IOobject::NO_WRITE,
217  IOobject::NO_REGISTER
218  ),
219  mesh_,
221  )
222 {
223  read(dict);
224 
225  // If the flow is multiphase
226  if (!multiphaseThermo_.empty())
227  {
228  Cps_.setSize(multiphaseThermo_.size());
229  kappas_.setSize(Cps_.size());
230  phaseNames_.setSize(Cps_.size());
231 
232  label phasei = 0;
233  forAllConstIters(multiphaseThermo_, iter)
234  {
235  const word& key = iter().keyword();
236 
237  const dictionary* subDictPtr = multiphaseThermo_.findDict(key);
238 
239  if (!subDictPtr)
240  {
242  << "Found non-dictionary entry " << iter()
243  << " in top-level dictionary " << multiphaseThermo_
244  << exit(FatalError);
245  }
246 
247  const dictionary& dict = *subDictPtr;
248 
249  phaseNames_[phasei] = key;
250 
251  Cps_.set
252  (
253  phasei,
255  (
256  "Cp",
258  dict
259  )
260  );
261 
262  kappas_.set
263  (
264  phasei,
265  new dimensionedScalar //[J/m/s/K]
266  (
267  "kappa",
269  dict
270  )
271  );
272 
273  ++phasei;
274  }
275 
276  phases_.setSize(phaseNames_.size());
277  forAll(phaseNames_, i)
278  {
279  phases_.set
280  (
281  i,
282  mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
283  );
284  }
285 
286  rhoCp_ = rho()*Cp();
287  rhoCp_.oldTime();
288  }
289  else
290  {
291  if (Cp_.value() == 0.0 || kappa_.value() == 0.0)
292  {
294  << " Multiphase thermo dictionary not found and Cp/kappa "
295  << " for single phase are zero. Please entry either"
296  << exit(FatalError);
297  }
298 
299  }
300 
301  // Force creation of transported field so any BCs using it can
302  // look it up
303  volScalarField& s = transportedField();
304  s.correctBoundaryConditions();
305 }
306 
307 
308 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
311 {}
312 
313 
314 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
315 
317 {
319 
320  dict.readIfPresent("phi", phiName_);
321  dict.readIfPresent("rho", rhoName_);
322 
323  schemesField_ = dict.getOrDefault("schemesField", fieldName_);
324 
325  dict.readIfPresent("nCorr", nCorr_);
326 
327  if (dict.found("fvOptions"))
328  {
329  fvOptions_.reset(dict.subDict("fvOptions"));
330  }
331 
332  return true;
333 }
334 
335 
337 {
338  volScalarField& s = transportedField();
339 
340  Log << type() << " execute: " << s.name() << endl;
341 
342  const surfaceScalarField& phi =
343  mesh_.lookupObject<surfaceScalarField>(phiName_);
344 
345  // Calculate the diffusivity
346  const volScalarField kappaEff("kappaEff", this->kappaEff());
347 
348  word divScheme("div(phi," + schemesField_ + ")");
349  word laplacianScheme
350  (
351  "laplacian(kappaEff," + schemesField_ + ")"
352  );
353 
354  // Set under-relaxation coeff
355  scalar relaxCoeff = 0;
356  mesh_.relaxEquation(schemesField_, relaxCoeff);
357 
358  if (phi.dimensions() == dimMass/dimTime)
359  {
360  rhoCp_ = rho()*Cp();
362 
363  for (label i = 0; i <= nCorr_; i++)
364  {
365  fvScalarMatrix sEqn
366  (
367  fvm::ddt(rhoCp_, s)
368  + fvm::div(rhoCpPhi, s, divScheme)
369  - fvm::Sp(fvc::ddt(rhoCp_) + fvc::div(rhoCpPhi), s)
370  - fvm::laplacian(kappaEff, s, laplacianScheme)
371  ==
372  fvOptions_(rhoCp_, s)
373  );
374 
375  sEqn.relax(relaxCoeff);
376 
377  fvOptions_.constrain(sEqn);
378 
379  sEqn.solve(schemesField_);
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 (label i = 0; i <= nCorr_; i++)
397  {
398  fvScalarMatrix sEqn
399  (
400  fvm::ddt(rhoCp, s)
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  sEqn.solve(schemesField_);
412  }
413  }
414  else
415  {
417  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
418  << "Dimensions should be " << dimMass/dimTime << " or "
420  }
422  Log << endl;
423 
424  return true;
425 }
426 
427 
429 {
430  return true;
431 }
432 
433 
434 // ************************************************************************* //
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:608
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
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:72
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:320
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
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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:1095
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
#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:180
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 it is a dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127