scalarTransport.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 | Copyright (C) 2012-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "scalarTransport.H"
27 #include "surfaceFields.H"
28 #include "dictionary.H"
31 #include "fvScalarMatrix.H"
32 #include "fvmDdt.H"
33 #include "fvmDiv.H"
34 #include "fvcDiv.H"
35 #include "fvmLaplacian.H"
36 #include "fvmSup.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(scalarTransport, 0);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 {
53 
54  wordList bTypes(U.boundaryField().size());
55 
56  forAll(bTypes, patchI)
57  {
58  const fvPatchField<vector>& pf = U.boundaryField()[patchI];
59  if (isA<fixedValueFvPatchVectorField>(pf))
60  {
61  bTypes[patchI] = fixedValueFvPatchScalarField::typeName;
62  }
63  else
64  {
65  bTypes[patchI] = zeroGradientFvPatchScalarField::typeName;
66  }
67  }
68 
69  return bTypes;
70 }
71 
72 
74 {
75  if (!mesh_.foundObject<volScalarField>(name()))
76  {
77  volScalarField* fldPtr = new volScalarField
78  (
79  IOobject
80  (
81  name(),
82  mesh_.time().timeName(),
83  mesh_,
86  ),
87  mesh_,
88  dimensionedScalar("zero", dimless, 0.0),
89  boundaryTypes()
90  );
91  fldPtr->store();
92  }
93 
94  return const_cast<volScalarField&>
95  (
96  mesh_.lookupObject<volScalarField>(name())
97  );
98 }
99 
100 
102 (
103  const surfaceScalarField& phi
104 ) const
105 {
106  typedef incompressible::turbulenceModel icoModel;
107  typedef compressible::turbulenceModel cmpModel;
108 
109  if (userDT_)
110  {
111  return tmp<volScalarField>
112  (
113  new volScalarField
114  (
115  IOobject
116  (
117  "DT",
118  mesh_.time().timeName(),
119  mesh_.time(),
122  ),
123  mesh_,
124  dimensionedScalar("DT", phi.dimensions()/dimLength, DT_)
125  )
126  );
127  }
128  else if (mesh_.foundObject<icoModel>(turbulenceModel::propertiesName))
129  {
130  const icoModel& model = mesh_.lookupObject<icoModel>
131  (
133  );
134 
135  return model.nuEff();
136  }
137  else if (mesh_.foundObject<cmpModel>(turbulenceModel::propertiesName))
138  {
139  const cmpModel& model = mesh_.lookupObject<cmpModel>
140  (
142  );
143 
144  return model.muEff();
145  }
146  else
147  {
148  return tmp<volScalarField>
149  (
150  new volScalarField
151  (
152  IOobject
153  (
154  "DT",
155  mesh_.time().timeName(),
156  mesh_.time(),
159  ),
160  mesh_,
161  dimensionedScalar("DT", phi.dimensions()/dimLength, 0.0)
162  )
163  );
164  }
165 }
166 
167 
168 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
169 
171 (
172  const word& name,
173  const objectRegistry& obr,
174  const dictionary& dict,
175  const bool loadFromFiles
176 )
177 :
178  name_(name),
179  mesh_(refCast<const fvMesh>(obr)),
180  active_(true),
181  phiName_("phi"),
182  UName_("U"),
183  rhoName_("rho"),
184  DT_(0.0),
185  userDT_(false),
186  resetOnStartUp_(false),
187  nCorr_(0),
188  autoSchemes_(false),
189  fvOptions_(mesh_),
190  log_(true)
191 {
192  read(dict);
193 
194  // Force creation of transported field so any bcs using it can look it
195  // up
196  volScalarField& T = transportedField();
197 
198  if (resetOnStartUp_)
199  {
200  T == dimensionedScalar("zero", dimless, 0.0);
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
206 
208 {}
209 
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
214 {
215  if (active_)
216  {
217  log_.readIfPresent("log", dict);
218 
219  if (log_) Info<< type() << " " << name_ << " output:" << nl;
220 
221  dict.readIfPresent("phiName", phiName_);
222  dict.readIfPresent("UName", UName_);
223  dict.readIfPresent("rhoName", rhoName_);
224 
225  userDT_ = false;
226  if (dict.readIfPresent("DT", DT_))
227  {
228  userDT_ = true;
229  }
230 
231  dict.readIfPresent("nCorr", nCorr_);
232  dict.lookup("resetOnStartUp") >> resetOnStartUp_;
233  dict.lookup("autoSchemes") >> autoSchemes_;
234 
235  fvOptions_.reset(dict.subDict("fvOptions"));
236  }
237 }
238 
239 
241 {
242  if (active_)
243  {
244  if (log_) Info<< type() << " " << name_ << " output:" << nl;
245 
246  const surfaceScalarField& phi =
247  mesh_.lookupObject<surfaceScalarField>(phiName_);
248 
249  volScalarField& T = transportedField();
250 
251  // calculate the diffusivity
252  volScalarField DT(this->DT(phi));
253 
254  // set schemes
255  word schemeVar = T.name();
256  if (autoSchemes_)
257  {
258  schemeVar = UName_;
259  }
260 
261  word divScheme("div(phi," + schemeVar + ")");
262  word laplacianScheme("laplacian(" + DT.name() + "," + schemeVar + ")");
263 
264  // set under-relaxation coeff
265  scalar relaxCoeff = 0.0;
266  if (mesh_.relaxEquation(schemeVar))
267  {
268  relaxCoeff = mesh_.equationRelaxationFactor(schemeVar);
269  }
270 
271  if (phi.dimensions() == dimMass/dimTime)
272  {
273  const volScalarField& rho =
274  mesh_.lookupObject<volScalarField>(rhoName_);
275 
276  // solve
277  for (label i = 0; i <= nCorr_; i++)
278  {
280  (
281  fvm::ddt(rho, T)
282  + fvm::div(phi, T, divScheme)
283  - fvm::laplacian(DT, T, laplacianScheme)
284  ==
285  fvOptions_(rho, T)
286  );
287 
288  TEqn.relax(relaxCoeff);
289 
290  fvOptions_.constrain(TEqn);
291 
292  TEqn.solve(mesh_.solverDict(schemeVar));
293  }
294  }
295  else if (phi.dimensions() == dimVolume/dimTime)
296  {
297  // solve
298  for (label i = 0; i <= nCorr_; i++)
299  {
301  (
302  fvm::ddt(T)
303  + fvm::div(phi, T, divScheme)
304  - fvm::laplacian(DT, T, laplacianScheme)
305  ==
306  fvOptions_(T)
307  );
308 
309  TEqn.relax(relaxCoeff);
310 
311  fvOptions_.constrain(TEqn);
312 
313  TEqn.solve(mesh_.solverDict(schemeVar));
314  }
315  }
316  else
317  {
319  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
320  << "Dimensions should be " << dimMass/dimTime << " or "
321  << dimVolume/dimTime << endl;
322  }
323 
324  if (log_) Info<< endl;
325  }
326 }
327 
328 
330 {
331  // Do nothing
332 }
333 
334 
336 {
337  // Do nothing
338 }
339 
340 
342 {
343  // Do nothing
344 }
345 
346 
347 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::scalarTransport::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: scalarTransport.H:175
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::scalarTransport::write
virtual void write()
Calculate the scalarTransport and write.
Definition: scalarTransport.C:341
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
turbulentTransportModel.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
fvcDiv.H
Calculate the divergence of the given field.
Foam::scalarTransport::boundaryTypes
wordList boundaryTypes() const
Return the boundary types for the scalar field.
Definition: scalarTransport.C:50
Foam::scalarTransport::UName_
word UName_
Name of velocity field (optional)
Definition: scalarTransport.H:184
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
U
U
Definition: pEqn.H:46
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Foam::scalarTransport::scalarTransport
scalarTransport(const scalarTransport &)
Disallow default bitwise copy construct.
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::scalarTransport::end
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: scalarTransport.C:329
Foam::scalarTransport::execute
virtual void execute()
Execute, currently does nothing.
Definition: scalarTransport.C:240
Foam::scalarTransport::transportedField
volScalarField & transportedField()
Return reference to registered transported field.
Definition: scalarTransport.C:73
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:519
fvmLaplacian.H
Calculate the matrix for the laplacian of the field.
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
fvScalarMatrix.H
A scalar instance of fvMatrix.
Foam::scalarTransport::DT
tmp< volScalarField > DT(const surfaceScalarField &phi) const
Return the diffusivity field.
Definition: scalarTransport.C:102
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
fvmSup.H
Calculate the matrix for implicit and explicit sources.
TEqn
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
rho
rho
Definition: pEqn.H:3
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
fixedValueFvPatchFields.H
dictionary.H
Foam::scalarTransport::read
virtual void read(const dictionary &)
Read the scalarTransport data.
Definition: scalarTransport.C:213
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::scalarTransport::timeSet
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: scalarTransport.C:335
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:52
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
scalarTransport.H
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
zeroGradientFvPatchFields.H
Foam::scalarTransport::~scalarTransport
virtual ~scalarTransport()
Destructor.
Definition: scalarTransport.C:207
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
turbulentFluidThermoModel.H
fvmDdt.H
Calulate the matrix for the first temporal derivative.