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-2020 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 {
62  typedef compressible::momentumTransportModel cmpModel;
63 
64  word Dname("D" + s_.name());
65 
66  if (constantD_)
67  {
68  return volScalarField::New
69  (
70  Dname,
71  mesh_,
72  dimensionedScalar(Dname, phi.dimensions()/dimLength, D_)
73  );
74  }
75  else if (mesh_.foundObject<icoModel>(momentumTransportModel::typeName))
76  {
77  const icoModel& model = mesh_.lookupObject<icoModel>
78  (
79  momentumTransportModel::typeName
80  );
81 
82  return alphaD_*model.nu() + alphaDt_*model.nut();
83  }
84  else if (mesh_.foundObject<cmpModel>(momentumTransportModel::typeName))
85  {
86  const cmpModel& model = mesh_.lookupObject<cmpModel>
87  (
88  momentumTransportModel::typeName
89  );
90 
91  return alphaD_*model.mu() + alphaDt_*model.mut();
92  }
93  else
94  {
95  return volScalarField::New
96  (
97  Dname,
98  mesh_,
99  dimensionedScalar(Dname, phi.dimensions()/dimLength, 0)
100  );
101  }
102 }
103 
104 
105 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
106 
108 (
109  const word& name,
110  const Time& runTime,
111  const dictionary& dict
112 )
113 :
114  fvMeshFunctionObject(name, runTime, dict),
115  fieldName_(dict.lookupOrDefault<word>("field", "s")),
116  D_(0),
117  nCorr_(0),
118  fvOptions_(mesh_),
119  s_
120  (
121  IOobject
122  (
123  fieldName_,
124  mesh_.time().timeName(),
125  mesh_,
128  ),
129  mesh_
130  )
131 {
132  read(dict);
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
139 {}
140 
141 
142 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
143 
145 {
147 
148  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
149  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
150  schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
151 
152  constantD_ = dict.readIfPresent("D", D_);
153  alphaD_ = dict.lookupOrDefault("alphaD", 1.0);
154  alphaDt_ = dict.lookupOrDefault("alphaDt", 1.0);
155 
156  dict.readIfPresent("nCorr", nCorr_);
157 
158  if (dict.found("fvOptions"))
159  {
160  fvOptions_.reset(dict.subDict("fvOptions"));
161  }
162 
163  return true;
164 }
165 
166 
168 {
169  Info<< type() << " write:" << endl;
170 
171  const surfaceScalarField& phi =
172  mesh_.lookupObject<surfaceScalarField>(phiName_);
173 
174  // Calculate the diffusivity
175  volScalarField D("D" + s_.name(), this->D(phi));
176 
177  word divScheme("div(phi," + schemesField_ + ")");
178  word laplacianScheme("laplacian(" + D.name() + "," + schemesField_ + ")");
179 
180  // Set under-relaxation coeff
181  scalar relaxCoeff = 0.0;
182  if (mesh_.relaxEquation(schemesField_))
183  {
184  relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
185  }
186 
187  if (phi.dimensions() == dimMass/dimTime)
188  {
189  const volScalarField& rho =
190  mesh_.lookupObject<volScalarField>(rhoName_);
191 
192  for (int i=0; i<=nCorr_; i++)
193  {
194  fvScalarMatrix sEqn
195  (
196  fvm::ddt(rho, s_)
197  + fvm::div(phi, s_, divScheme)
198  - fvm::laplacian(D, s_, laplacianScheme)
199  ==
200  fvOptions_(rho, s_)
201  );
202 
203  sEqn.relax(relaxCoeff);
204 
205  fvOptions_.constrain(sEqn);
206 
207  sEqn.solve(schemesField_);
208  }
209  }
210  else if (phi.dimensions() == dimVolume/dimTime)
211  {
212  for (int i=0; i<=nCorr_; i++)
213  {
214  fvScalarMatrix sEqn
215  (
216  fvm::ddt(s_)
217  + fvm::div(phi, s_, divScheme)
218  - fvm::laplacian(D, s_, laplacianScheme)
219  ==
220  fvOptions_(s_)
221  );
222 
223  sEqn.relax(relaxCoeff);
224 
225  fvOptions_.constrain(sEqn);
226 
227  sEqn.solve(schemesField_);
228  }
229  }
230  else
231  {
233  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
234  << "Dimensions should be " << dimMass/dimTime << " or "
236  }
237 
238  Info<< endl;
239 
240  return true;
241 }
242 
243 
245 {
246  return true;
247 }
248 
249 
250 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
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:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
IncompressibleMomentumTransportModel< transportModel > momentumTransportModel
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Calculate the matrix for the laplacian of the field.
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.
scalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
virtual bool execute()
Calculate the scalarTransport.
CompressibleMomentumTransportModel< fluidThermo > momentumTransportModel
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
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:58
static const char nl
Definition: Ostream.H:260
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:521
defineTypeNameAndDebug(Qdot, 0)
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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...
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.
Namespace for OpenFOAM.