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-2016 OpenFOAM Foundation
6  \\/ M anipulation |
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 
50 Foam::wordList Foam::scalarTransport::boundaryTypes() const
51 {
52  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
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 
73 Foam::tmp<Foam::volScalarField> Foam::scalarTransport::DT
74 (
75  const surfaceScalarField& phi
76 ) const
77 {
78  typedef incompressible::turbulenceModel icoModel;
79  typedef compressible::turbulenceModel cmpModel;
80 
81  if (userDT_)
82  {
83  return tmp<volScalarField>
84  (
85  new volScalarField
86  (
87  IOobject
88  (
89  "DT",
90  mesh_.time().timeName(),
91  mesh_.time(),
94  ),
95  mesh_,
96  dimensionedScalar("DT", phi.dimensions()/dimLength, DT_)
97  )
98  );
99  }
100  else if (mesh_.foundObject<icoModel>(turbulenceModel::propertiesName))
101  {
102  const icoModel& model = mesh_.lookupObject<icoModel>
103  (
105  );
106 
107  return model.nuEff();
108  }
109  else if (mesh_.foundObject<cmpModel>(turbulenceModel::propertiesName))
110  {
111  const cmpModel& model = mesh_.lookupObject<cmpModel>
112  (
114  );
115 
116  return model.muEff();
117  }
118  else
119  {
120  return tmp<volScalarField>
121  (
122  new volScalarField
123  (
124  IOobject
125  (
126  "DT",
127  mesh_.time().timeName(),
128  mesh_.time(),
131  ),
132  mesh_,
133  dimensionedScalar("DT", phi.dimensions()/dimLength, 0.0)
134  )
135  );
136  }
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
142 Foam::scalarTransport::scalarTransport
143 (
144  const word& name,
145  const objectRegistry& obr,
146  const dictionary& dict,
147  const bool loadFromFiles
148 )
149 :
150  name_(name),
151  mesh_(refCast<const fvMesh>(obr)),
152  phiName_(dict.lookupOrDefault<word>("phiName", "phi")),
153  UName_(dict.lookupOrDefault<word>("UName", "U")),
154  rhoName_(dict.lookupOrDefault<word>("rhoName", "rho")),
155  DT_(0.0),
156  userDT_(false),
157  resetOnStartUp_(false),
158  nCorr_(0),
159  autoSchemes_(false),
160  fvOptions_(mesh_),
161  T_
162  (
163  IOobject
164  (
165  name,
166  mesh_.time().timeName(),
167  mesh_,
170  ),
171  mesh_,
172  dimensionedScalar("zero", dimless, 0.0),
173  boundaryTypes()
174  )
175 {
176  read(dict);
177 
178  if (resetOnStartUp_)
179  {
180  T_ == dimensionedScalar("zero", dimless, 0.0);
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
186 
188 {}
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
194 {
195  Info<< type() << ":" << nl;
196 
197  phiName_ = dict.lookupOrDefault<word>("phiName", "phi");
198  UName_ = dict.lookupOrDefault<word>("UName", "U");
199  rhoName_ = dict.lookupOrDefault<word>("rhoName", "rho");
200 
201  userDT_ = false;
202  if (dict.readIfPresent("DT", DT_))
203  {
204  userDT_ = true;
205  }
206 
207  dict.lookup("resetOnStartUp") >> resetOnStartUp_;
208 
209  dict.readIfPresent("nCorr", nCorr_);
210 
211  dict.lookup("autoSchemes") >> autoSchemes_;
212 
213  fvOptions_.reset(dict.subDict("fvOptions"));
214 }
215 
216 
218 {
219  Info<< type() << " output:" << endl;
220 
221  const surfaceScalarField& phi =
222  mesh_.lookupObject<surfaceScalarField>(phiName_);
223 
224  // calculate the diffusivity
225  volScalarField DT(this->DT(phi));
226 
227  // set schemes
228  word schemeVar = T_.name();
229  if (autoSchemes_)
230  {
231  schemeVar = UName_;
232  }
233 
234  word divScheme("div(phi," + schemeVar + ")");
235  word laplacianScheme("laplacian(" + DT.name() + "," + schemeVar + ")");
236 
237  // set under-relaxation coeff
238  scalar relaxCoeff = 0.0;
239  if (mesh_.relaxEquation(schemeVar))
240  {
241  relaxCoeff = mesh_.equationRelaxationFactor(schemeVar);
242  }
243 
244  if (phi.dimensions() == dimMass/dimTime)
245  {
246  const volScalarField& rho =
247  mesh_.lookupObject<volScalarField>(rhoName_);
248 
249  // solve
250  for (label i = 0; i <= nCorr_; i++)
251  {
253  (
254  fvm::ddt(rho, T_)
255  + fvm::div(phi, T_, divScheme)
256  - fvm::laplacian(DT, T_, laplacianScheme)
257  ==
258  fvOptions_(rho, T_)
259  );
260 
261  TEqn.relax(relaxCoeff);
262 
263  fvOptions_.constrain(TEqn);
264 
265  TEqn.solve(mesh_.solverDict(schemeVar));
266  }
267  }
268  else if (phi.dimensions() == dimVolume/dimTime)
269  {
270  // solve
271  for (label i = 0; i <= nCorr_; i++)
272  {
274  (
275  fvm::ddt(T_)
276  + fvm::div(phi, T_, divScheme)
277  - fvm::laplacian(DT, T_, laplacianScheme)
278  ==
279  fvOptions_(T_)
280  );
281 
282  TEqn.relax(relaxCoeff);
283 
284  fvOptions_.constrain(TEqn);
285 
286  TEqn.solve(mesh_.solverDict(schemeVar));
287  }
288  }
289  else
290  {
292  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
293  << "Dimensions should be " << dimMass/dimTime << " or "
294  << dimVolume/dimTime << endl;
295  }
296 
297  Info<< endl;
298 }
299 
300 
302 {
303  execute();
304 }
305 
306 
308 {}
309 
310 
312 {}
313 
314 
315 // ************************************************************************* //
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::surfaceFields.
virtual ~scalarTransport()
Destructor.
Calculate the matrix for implicit and explicit sources.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T)-fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for handling words, derived from string.
Definition: word.H:59
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
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Calculate the matrix for the divergence of the given field and flux.
virtual void execute()
Execute, currently does nothing.
virtual void write()
Calculate the scalarTransport and write.
Calulate the matrix for the first temporal derivative.
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:68
messageStream Info
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Namespace for OpenFOAM.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
#define forAll(list, i)
Definition: UList.H:421
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:525
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const dimensionSet & dimensions() const
Return dimensions.
solverPerformance solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:57
virtual void end()
Execute at the final time-loop, currently does nothing.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Calculate the matrix for the laplacian of the field.
Calculate the divergence of the given field.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
List< word > wordList
A List of words.
Definition: fileName.H:54
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Registry of regIOobjects.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual void read(const dictionary &)
Read the scalarTransport data.
static const word propertiesName
Default name of the turbulence properties dictionary.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
A scalar instance of fvMatrix.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
IncompressibleTurbulenceModel< transportModel > turbulenceModel
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
defineTypeNameAndDebug(combustionModel, 0)