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-2020 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  ),
125  mesh_,
126  rho_
127  );
128 
129  if (phases_.size())
130  {
131  trho.ref() = lookupObject<volScalarField>(rhoName_);
132  }
133  return trho;
134 }
135 
136 
138 Foam::functionObjects::energyTransport::Cp() const
139 {
140  if (phases_.size())
141  {
142  tmp<volScalarField> tCp(phases_[0]*Cps_[0]);
143 
144  for (label i = 1; i < phases_.size(); i++)
145  {
146  tCp.ref() += phases_[i]*Cps_[i];
147  }
148  return tCp;
149  }
150 
152  (
153  IOobject
154  (
155  "tCp",
156  mesh_.time().timeName(),
157  mesh_,
160  ),
161  mesh_,
162  Cp_
163  );
164 }
165 
166 
168 Foam::functionObjects::energyTransport::kappa() const
169 {
170  if (phases_.size())
171  {
172  tmp<volScalarField> tkappa(phases_[0]*kappas_[0]);
173 
174  for (label i = 1; i < phases_.size(); i++)
175  {
176  tkappa.ref() += phases_[i]*kappas_[i];
177  }
178  return tkappa;
179  }
180 
182  (
183  IOobject
184  (
185  "tkappa",
186  mesh_.time().timeName(),
187  mesh_,
190  ),
191  mesh_,
192  kappa_
193  );
194 }
195 
196 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
197 
198 Foam::functionObjects::energyTransport::energyTransport
199 (
200  const word& name,
201  const Time& runTime,
202  const dictionary& dict
203 )
204 :
206  fieldName_(dict.getOrDefault<word>("field", "T")),
207  phiName_(dict.getOrDefault<word>("phi", "phi")),
208  rhoName_(dict.getOrDefault<word>("rho", "rho")),
209  nCorr_(0),
210  schemesField_("unknown-schemesField"),
211  fvOptions_(mesh_),
212  multiphaseThermo_(dict.subOrEmptyDict("phaseThermos")),
213  Cp_("Cp", dimEnergy/dimMass/dimTemperature, 0, dict),
214  kappa_
215  (
216  "kappa",
218  0,
219  dict
220  ),
221  rho_("rhoInf", dimDensity, 0, dict),
222  Prt_("Prt", dimless, 1, dict),
223  rhoCp_
224  (
225  IOobject
226  (
227  "rhoCp",
228  mesh_.time().timeName(),
229  mesh_,
232  ),
233  mesh_,
235  )
236 {
237  read(dict);
238 
239  // If the flow is multiphase
240  if (!multiphaseThermo_.empty())
241  {
242  Cps_.setSize(multiphaseThermo_.size());
243  kappas_.setSize(Cps_.size());
244  phaseNames_.setSize(Cps_.size());
245 
246  label phasei = 0;
247  forAllConstIters(multiphaseThermo_, iter)
248  {
249  const word& key = iter().keyword();
250 
251  if (!multiphaseThermo_.isDict(key))
252  {
254  << "Found non-dictionary entry " << iter()
255  << " in top-level dictionary " << multiphaseThermo_
256  << exit(FatalError);
257  }
258 
259  const dictionary& dict = multiphaseThermo_.subDict(key);
260 
261  phaseNames_[phasei] = key;
262 
263  Cps_.set
264  (
265  phasei,
267  (
268  "Cp",
270  dict
271  )
272  );
273 
274  kappas_.set
275  (
276  phasei,
277  new dimensionedScalar //[J/m/s/K]
278  (
279  "kappa",
281  dict
282  )
283  );
284 
285  ++phasei;
286  }
287 
288  phases_.setSize(phaseNames_.size());
289  forAll(phaseNames_, i)
290  {
291  phases_.set
292  (
293  i,
294  mesh_.getObjectPtr<volScalarField>(phaseNames_[i])
295  );
296  }
297 
298  rhoCp_ = rho()*Cp();
299  rhoCp_.oldTime();
300  }
301  else
302  {
303  if (Cp_.value() == 0.0 || kappa_.value() == 0.0)
304  {
306  << " Multiphase thermo dictionary not found and Cp/kappa "
307  << " for single phase are zero. Please entry either"
308  << exit(FatalError);
309  }
310 
311  }
312 
313  // Force creation of transported field so any BCs using it can
314  // look it up
315  volScalarField& s = transportedField();
316  s.correctBoundaryConditions();
317 }
318 
319 
320 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
321 
323 {}
324 
325 
326 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
327 
329 {
331 
332  dict.readIfPresent("phi", phiName_);
333  dict.readIfPresent("rho", rhoName_);
334 
335  schemesField_ = dict.getOrDefault("schemesField", fieldName_);
336 
337  dict.readIfPresent("nCorr", nCorr_);
338 
339  if (dict.found("fvOptions"))
340  {
341  fvOptions_.reset(dict.subDict("fvOptions"));
342  }
343 
344  return true;
345 }
346 
347 
349 {
350  volScalarField& s = transportedField();
351 
352  Log << type() << " execute: " << s.name() << endl;
353 
354  const surfaceScalarField& phi =
355  mesh_.lookupObject<surfaceScalarField>(phiName_);
356 
357  // Calculate the diffusivity
358  const volScalarField kappaEff("kappaEff", this->kappaEff());
359 
360  word divScheme("div(phi," + schemesField_ + ")");
361  word laplacianScheme
362  (
363  "laplacian(kappaEff," + schemesField_ + ")"
364  );
365 
366  // Set under-relaxation coeff
367  scalar relaxCoeff = 0.0;
368  if (mesh_.relaxEquation(schemesField_))
369  {
370  relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
371  }
372 
373  if (phi.dimensions() == dimMass/dimTime)
374  {
375  rhoCp_ = rho()*Cp();
377 
378  for (label i = 0; i <= nCorr_; i++)
379  {
380  fvScalarMatrix sEqn
381  (
382  fvm::ddt(rhoCp_, s)
383  + fvm::div(rhoCpPhi, s, divScheme)
384  - fvm::Sp(fvc::ddt(rhoCp_) + fvc::div(rhoCpPhi), s)
385  - fvm::laplacian(kappaEff, s, laplacianScheme)
386  ==
387  fvOptions_(rhoCp_, s)
388  );
389 
390  sEqn.relax(relaxCoeff);
391 
392  fvOptions_.constrain(sEqn);
393 
394  sEqn.solve(mesh_.solverDict(schemesField_));
395  }
396  }
397  else if (phi.dimensions() == dimVolume/dimTime)
398  {
399  dimensionedScalar rhoCp(rho_*Cp_);
400 
401  const surfaceScalarField CpPhi(rhoCp*phi);
402 
403  auto trhoCp = tmp<volScalarField>::New
404  (
405  IOobject
406  (
407  "trhoCp",
408  mesh_.time().timeName(),
409  mesh_,
412  ),
413  mesh_,
414  rhoCp
415  );
416 
417  for (label i = 0; i <= nCorr_; i++)
418  {
419  fvScalarMatrix sEqn
420  (
421  fvm::ddt(rhoCp, s)
422  + fvm::div(CpPhi, s, divScheme)
423  - fvm::laplacian(kappaEff, s, laplacianScheme)
424  ==
425  fvOptions_(trhoCp.ref(), s)
426  );
427 
428  sEqn.relax(relaxCoeff);
429 
430  fvOptions_.constrain(sEqn);
431 
432  sEqn.solve(mesh_.solverDict(schemesField_));
433  }
434  }
435  else
436  {
438  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
439  << "Dimensions should be " << dimMass/dimTime << " or "
441  }
442 
443  Log << endl;
444 
445  return true;
446 }
447 
448 
450 {
451  return true;
452 }
453 
454 
455 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Log
#define Log
Definition: PDRblock.C:35
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
s
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))
Definition: gmvOutputSpray.H:25
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
Foam::functionObjects::energyTransport::read
virtual bool read(const dictionary &)
Read the energyTransport data.
Definition: energyTransport.C:328
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
rhoCp
rhoCp
Definition: TEqn.H:3
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
kappaEff
kappaEff
Definition: TEqn.H:10
turbulentTransportModel.H
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
tCp
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
phasei
label phasei
Definition: pEqn.H:27
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
surfaceFields.H
Foam::surfaceFields.
Foam::functionObjects::energyTransport::execute
virtual bool execute()
Calculate the energyTransport.
Definition: energyTransport.C:348
Foam::functionObjects::regionFunctionObject::store
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
Definition: regionFunctionObjectTemplates.C:107
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::functionObjects::energyTransport::write
virtual bool write()
Do nothing.
Definition: energyTransport.C:449
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
rho
rho
Definition: readInitialConditions.H:88
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:559
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::fvc::interpolate
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.
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::nl
constexpr char nl
Definition: Ostream.H:385
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
Foam::functionObjects::energyTransport::~energyTransport
virtual ~energyTransport()
Destructor.
Definition: energyTransport.C:322
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:301
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:55
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
energyTransport.H
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
Foam::incompressible::turbulenceModel
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Definition: turbulentTransportModel.H:60
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
rhoCpPhi
const surfaceScalarField rhoCpPhi(fvc::interpolate(fluid.Cp()) *rhoPhi)
turbulentFluidThermoModel.H
fvmDdt.H
Calculate the matrix for the first temporal derivative.
Foam::IOobject::MUST_READ
Definition: IOobject.H:120