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-2021 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"
32 #include "fvModels.H"
33 #include "fvConstraints.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace functionObjects
43 {
44  defineTypeNameAndDebug(scalarTransport, 0);
45 
47  (
48  functionObject,
49  scalarTransport,
50  dictionary
51  );
52 }
53 }
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
58 Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
59 (
60  const surfaceScalarField& phi
61 ) const
62 {
64  typedef compressible::momentumTransportModel cmpModel;
65 
66  word Dname("D" + s_.name());
67 
68  if (constantD_)
69  {
70  return volScalarField::New
71  (
72  Dname,
73  mesh_,
74  dimensionedScalar(Dname, phi.dimensions()/dimLength, D_)
75  );
76  }
77  else if (mesh_.foundObject<icoModel>(momentumTransportModel::typeName))
78  {
79  const icoModel& model = mesh_.lookupObject<icoModel>
80  (
81  momentumTransportModel::typeName
82  );
83 
84  return alphaD_*model.nu() + alphaDt_*model.nut();
85  }
86  else if (mesh_.foundObject<cmpModel>(momentumTransportModel::typeName))
87  {
88  const cmpModel& model = mesh_.lookupObject<cmpModel>
89  (
90  momentumTransportModel::typeName
91  );
92 
93  return alphaD_*model.mu() + alphaDt_*model.mut();
94  }
95  else
96  {
97  return volScalarField::New
98  (
99  Dname,
100  mesh_,
101  dimensionedScalar(Dname, phi.dimensions()/dimLength, 0)
102  );
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
110 (
111  const word& name,
112  const Time& runTime,
113  const dictionary& dict
114 )
115 :
116  fvMeshFunctionObject(name, runTime, dict),
117  fieldName_(dict.lookupOrDefault<word>("field", "s")),
118  D_(0),
119  nCorr_(0),
120  s_
121  (
122  IOobject
123  (
124  fieldName_,
125  mesh_.time().timeName(),
126  mesh_,
129  ),
130  mesh_
131  )
132 {
133  read(dict);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
138 
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
146 {
148 
149  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
150  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
151  schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
152 
153  constantD_ = dict.readIfPresent("D", D_);
154  alphaD_ = dict.lookupOrDefault("alphaD", 1.0);
155  alphaDt_ = dict.lookupOrDefault("alphaDt", 1.0);
156 
157  dict.readIfPresent("nCorr", nCorr_);
158 
159  return true;
160 }
161 
162 
164 {
165  Info<< type() << " write:" << endl;
166 
167  const surfaceScalarField& phi =
168  mesh_.lookupObject<surfaceScalarField>(phiName_);
169 
170  // Calculate the diffusivity
171  volScalarField D("D" + s_.name(), this->D(phi));
172 
173  word divScheme("div(phi," + schemesField_ + ")");
174  word laplacianScheme("laplacian(" + D.name() + "," + schemesField_ + ")");
175 
176  // Set under-relaxation coeff
177  scalar relaxCoeff = 0.0;
178  if (mesh_.relaxEquation(schemesField_))
179  {
180  relaxCoeff = mesh_.equationRelaxationFactor(schemesField_);
181  }
182 
185  (
187  );
188 
189  if (phi.dimensions() == dimMass/dimTime)
190  {
191  const volScalarField& rho =
192  mesh_.lookupObject<volScalarField>(rhoName_);
193 
194  for (int i=0; i<=nCorr_; i++)
195  {
196  fvScalarMatrix sEqn
197  (
198  fvm::ddt(rho, s_)
199  + fvm::div(phi, s_, divScheme)
200  - fvm::laplacian(D, s_, laplacianScheme)
201  ==
202  fvModels.source(rho, s_)
203  );
204 
205  sEqn.relax(relaxCoeff);
206 
207  fvConstraints.constrain(sEqn);
208 
209  sEqn.solve(schemesField_);
210 
211  fvConstraints.constrain(s_);
212  }
213  }
214  else if (phi.dimensions() == dimVolume/dimTime)
215  {
216  for (int i=0; i<=nCorr_; i++)
217  {
218  fvScalarMatrix sEqn
219  (
220  fvm::ddt(s_)
221  + fvm::div(phi, s_, divScheme)
222  - fvm::laplacian(D, s_, laplacianScheme)
223  ==
224  fvModels.source(s_)
225  );
226 
227  sEqn.relax(relaxCoeff);
228 
229  fvConstraints.constrain(sEqn);
230 
231  sEqn.solve(schemesField_);
232 
233  fvConstraints.constrain(s_);
234  }
235  }
236  else
237  {
239  << "Incompatible dimensions for phi: " << phi.dimensions() << nl
240  << "Dimensions should be " << dimMass/dimTime << " or "
242  }
243 
244  Info<< endl;
245 
246  return true;
247 }
248 
249 
251 {
252  return true;
253 }
254 
255 
256 // ************************************************************************* //
Foam::surfaceFields.
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
IncompressibleMomentumTransportModel< kinematicTransportModel > 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.
const dimensionSet dimLength
scalarTransport(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual bool execute()
Calculate the scalarTransport.
CompressibleMomentumTransportModel< dynamicTransportModel > momentumTransportModel
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
Finite volume models.
Definition: fvModels.H:60
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
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
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:521
const dimensionSet dimMass
defineTypeNameAndDebug(Qdot, 0)
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
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,.
Finite volume constraints.
Definition: fvConstraints.H:57
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 dimVolume
messageStream Info
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
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.