scalarTransport.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2018 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 "fvmDdt.H"
29 #include "fvmDiv.H"
30 #include "fvmLaplacian.H"
31 #include "fvmSup.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace functionObjects
41 {
42  defineTypeNameAndDebug(scalarTransport, 0);
43 
45  (
46  functionObject,
47  scalarTransport,
48  dictionary
49  );
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
57 (
58  const surfaceScalarField& phi
59 ) const
60 {
61  typedef incompressible::turbulenceModel icoModel;
62  typedef compressible::turbulenceModel cmpModel;
63 
64  word Dname("D" + s_.name());
65 
66  if (constantD_)
67  {
68  return tmp<volScalarField>
69  (
70  new volScalarField
71  (
72  IOobject
73  (
74  Dname,
75  mesh_.time().timeName(),
76  mesh_.time(),
79  ),
80  mesh_,
81  dimensionedScalar(Dname, phi.dimensions()/dimLength, D_)
82  )
83  );
84  }
85  else if (mesh_.foundObject<icoModel>(turbulenceModel::propertiesName))
86  {
87  const icoModel& model = mesh_.lookupObject<icoModel>
88  (
90  );
91 
92  return alphaD_*model.nu() + alphaDt_*model.nut();
93  }
94  else if (mesh_.foundObject<cmpModel>(turbulenceModel::propertiesName))
95  {
96  const cmpModel& model = mesh_.lookupObject<cmpModel>
97  (
99  );
100 
101  return alphaD_*model.mu() + alphaDt_*model.mut();
102  }
103  else
104  {
105  return tmp<volScalarField>
106  (
107  new volScalarField
108  (
109  IOobject
110  (
111  Dname,
112  mesh_.time().timeName(),
113  mesh_.time(),
116  ),
117  mesh_,
118  dimensionedScalar(Dname, phi.dimensions()/dimLength, 0.0)
119  )
120  );
121  }
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
127 Foam::functionObjects::scalarTransport::scalarTransport
128 (
129  const word& name,
130  const Time& runTime,
131  const dictionary& dict
132 )
133 :
134  fvMeshFunctionObject(name, runTime, dict),
135  fieldName_(dict.lookupOrDefault<word>("field", "s")),
136  D_(0),
137  nCorr_(0),
138  fvOptions_(mesh_),
139  s_
140  (
141  IOobject
142  (
143  fieldName_,
144  mesh_.time().timeName(),
145  mesh_,
148  ),
149  mesh_
150  )
151 {
152  read(dict);
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
157 
159 {}
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
165 {
167 
168  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
169  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
170  schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
171 
172  constantD_ = dict.readIfPresent("D", D_);
173  alphaD_ = dict.lookupOrDefault("alphaD", 1.0);
174  alphaDt_ = dict.lookupOrDefault("alphaDt", 1.0);
175 
176  dict.readIfPresent("nCorr", nCorr_);
177 
178  if (dict.found("fvOptions"))
179  {
180  fvOptions_.reset(dict.subDict("fvOptions"));
181  }
182 
183  return true;
184 }
185 
186 
188 {
189  Info<< type() << " write:" << endl;
190 
191  const surfaceScalarField& phi =
192  mesh_.lookupObject<surfaceScalarField>(phiName_);
193 
194  // Calculate the diffusivity
195  volScalarField D(this->D(phi));
196 
197  word divScheme("div(phi," + schemesField_ + ")");
198  word laplacianScheme("laplacian(" + D.name() + "," + schemesField_ + ")");
199 
200  // Set under-relaxation coeff
201  scalar relaxCoeff = 0.0;
202  if (mesh_.relaxEquation(schemesField_))
203  {
204  relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
205  }
206 
207  if (phi.dimensions() == dimMass/dimTime)
208  {
209  const volScalarField& rho =
210  mesh_.lookupObject<volScalarField>(rhoName_);
211 
212  for (label i = 0; i <= nCorr_; i++)
213  {
214  fvScalarMatrix sEqn
215  (
216  fvm::ddt(rho, s_)
217  + fvm::div(phi, s_, divScheme)
218  - fvm::laplacian(D, s_, laplacianScheme)
219  ==
220  fvOptions_(rho, s_)
221  );
222 
223  sEqn.relax(relaxCoeff);
224 
225  fvOptions_.constrain(sEqn);
226 
227  sEqn.solve(mesh_.solverDict(schemesField_));
228  }
229  }
230  else if (phi.dimensions() == dimVolume/dimTime)
231  {
232  for (label i = 0; i <= nCorr_; i++)
233  {
234  fvScalarMatrix sEqn
235  (
236  fvm::ddt(s_)
237  + fvm::div(phi, s_, divScheme)
238  - fvm::laplacian(D, s_, laplacianScheme)
239  ==
240  fvOptions_(s_)
241  );
242 
243  sEqn.relax(relaxCoeff);
244 
245  fvOptions_.constrain(sEqn);
246 
247  sEqn.solve(mesh_.solverDict(schemesField_));
248  }
249  }
250  else
251  {
253  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
254  << "Dimensions should be " << dimMass/dimTime << " or "
256  }
257 
258  Info<< endl;
259 
260  return true;
261 }
262 
263 
265 {
266  return true;
267 }
268 
269 
270 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Calculate the matrix for the laplacian of the field.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
virtual bool execute()
Calculate the scalarTransport.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual bool read(const dictionary &)
Read the scalarTransport data.
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:72
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:57
static const char nl
Definition: Ostream.H:265
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:523
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
Calculate the matrix for the divergence of the given field and flux.
virtual bool write()
Do nothing.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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:46
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Calculate the matrix for implicit and explicit sources.
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.