volumeFractionSource.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) 2019-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 "volumeFractionSource.H"
27 #include "fvmDiv.H"
28 #include "fvmLaplacian.H"
29 #include "fvcDiv.H"
30 #include "surfaceInterpolate.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace fv
39 {
40  defineTypeNameAndDebug(volumeFractionSource, 0);
42  (
43  fvModel,
44  volumeFractionSource,
45  dictionary
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::fv::volumeFractionSource::readCoeffs()
54 {
55  phiName_ = coeffs().lookupOrDefault<word>("phi", "phi");
56  rhoName_ = coeffs().lookupOrDefault<word>("rho", "rho");
57  UName_ = coeffs().lookupOrDefault<word>("U", "U");
58 
59  volumePhaseName_ = coeffs().lookup<word>("volumePhase");
60 }
61 
62 
63 const Foam::volScalarField& Foam::fv::volumeFractionSource::volumeAlpha() const
64 {
65  const word alphaName = IOobject::groupName("alpha", volumePhaseName_);
66 
67  if (!mesh().foundObject<volScalarField>(alphaName))
68  {
69  volScalarField* alphaPtr =
70  new volScalarField
71  (
72  IOobject
73  (
74  alphaName,
75  mesh().time().constant(),
76  mesh(),
79  ),
80  mesh()
81  );
82 
83  alphaPtr->store();
84  }
85 
86  return mesh().lookupObject<volScalarField>(alphaName);
87 }
88 
89 
90 Foam::tmp<Foam::volScalarField> Foam::fv::volumeFractionSource::D
91 (
92  const word& fieldName
93 ) const
94 {
95  const word phiName =
96  IOobject::groupName(phiName_, IOobject::group(fieldName));
97  const surfaceScalarField& phi =
98  mesh().lookupObject<surfaceScalarField>(phiName);
99 
100  if (phi.dimensions() == dimVolume/dimTime)
101  {
102  const momentumTransportModel& turbulence =
103  mesh().lookupObject<momentumTransportModel>
104  (
105  momentumTransportModel::typeName
106  );
107 
108  return turbulence.nuEff();
109  }
110  else if (phi.dimensions() == dimMass/dimTime)
111  {
112  const thermophysicalTransportModel& ttm =
113  mesh().lookupObject<thermophysicalTransportModel>
114  (
115  thermophysicalTransportModel::typeName
116  );
117 
118  return
119  fieldName == ttm.thermo().T().name()
120  ? ttm.kappaEff()
121  : fieldName == ttm.thermo().he().name()
122  ? ttm.alphaEff()
123  : ttm.momentumTransport().muEff();
124  }
125  else
126  {
128  << "Dimensions of " << phi.name() << " not recognised"
129  << exit(FatalError);
130  return tmp<volScalarField>(nullptr);
131  }
132 }
133 
134 
135 template <class Type, class AlphaFieldType>
136 void Foam::fv::volumeFractionSource::addGeneralSup
137 (
138  const AlphaFieldType& alpha,
139  fvMatrix<Type>& eqn,
140  const word& fieldName
141 ) const
142 {
143  const word phiName =
144  IOobject::groupName(phiName_, IOobject::group(fieldName));
145  const surfaceScalarField& phi =
146  mesh().lookupObject<surfaceScalarField>(phiName);
147 
148  const volScalarField B(1 - volumeAlpha());
149  const volScalarField AByB(volumeAlpha()/B);
150  const volScalarField D(this->D(fieldName));
151 
152  // Divergence term
153  const word divScheme = "div(" + phiName + "," + eqn.psi().name() + ")";
154  eqn -= AByB*fvm::div(phi, eqn.psi(), divScheme);
155 
156  // Laplacian term
157  const word laplacianScheme =
158  "laplacian(" + D.name() + "," + eqn.psi().name() + ")";
159  eqn +=
160  fvm::laplacian(D, eqn.psi())
161  - 1/B*fvm::laplacian(B*D, eqn.psi(), laplacianScheme);
162 }
163 
164 
165 template<class Type, class AlphaFieldType>
166 void Foam::fv::volumeFractionSource::addAlphaSupType
167 (
168  const AlphaFieldType& alpha,
169  fvMatrix<Type>& eqn,
170  const word& fieldName
171 ) const
172 {
173  addGeneralSup(alpha, eqn, fieldName);
174 }
175 
176 
177 template<class AlphaFieldType>
178 void Foam::fv::volumeFractionSource::addAlphaSupType
179 (
180  const AlphaFieldType& alpha,
181  fvMatrix<scalar>& eqn,
182  const word& fieldName
183 ) const
184 {
185  if (IOobject::member(fieldName) == rhoName_)
186  {
187  const word phiName =
188  IOobject::groupName(phiName_, IOobject::group(fieldName));
189  const surfaceScalarField& phi =
190  mesh().lookupObject<surfaceScalarField>(phiName);
191 
192  const volScalarField AByB(volumeAlpha()/(1 - volumeAlpha()));
193 
194  eqn -= AByB*fvc::div(phi);
195  }
196  else
197  {
198  addGeneralSup(alpha, eqn, fieldName);
199  }
200 }
201 
202 
203 template<class AlphaFieldType>
204 void Foam::fv::volumeFractionSource::addAlphaSupType
205 (
206  const AlphaFieldType& alpha,
207  fvMatrix<vector>& eqn,
208  const word& fieldName
209 ) const
210 {
211  if (IOobject::member(fieldName) == UName_)
212  {
213  const word phiName =
214  IOobject::groupName(phiName_, IOobject::group(fieldName));
215  const surfaceScalarField& phi =
216  mesh().lookupObject<surfaceScalarField>(phiName);
217 
218  const volScalarField AByB(volumeAlpha()/(1 - volumeAlpha()));
219 
220  const word scheme("div(" + phiName + "," + eqn.psi().name() + ")");
221 
222  eqn -= fvm::div(fvc::interpolate(AByB)*phi, eqn.psi(), scheme);
223  }
224  else
225  {
226  addGeneralSup(alpha, eqn, fieldName);
227  }
228 }
229 
230 
231 template<class Type>
232 void Foam::fv::volumeFractionSource::addSupType
233 (
234  fvMatrix<Type>& eqn,
235  const word& fieldName
236 ) const
237 {
238  addAlphaSupType(geometricOneField(), eqn, fieldName);
239 }
240 
241 
242 template<class Type>
243 void Foam::fv::volumeFractionSource::addSupType
244 (
245  const volScalarField& rho,
246  fvMatrix<Type>& eqn,
247  const word& fieldName
248 ) const
249 {
250  addAlphaSupType(geometricOneField(), eqn, fieldName);
251 }
252 
253 
254 template<class Type>
255 void Foam::fv::volumeFractionSource::addSupType
256 (
257  const volScalarField& alpha,
258  const volScalarField& rho,
259  fvMatrix<Type>& eqn,
260  const word& fieldName
261 ) const
262 {
263  addAlphaSupType(alpha, eqn, fieldName);
264 }
265 
266 
267 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
268 
270 (
271  const word& name,
272  const word& modelType,
273  const dictionary& dict,
274  const fvMesh& mesh
275 )
276 :
277  fvModel(name, modelType, dict, mesh),
278  phiName_(word::null),
279  rhoName_(word::null),
280  UName_(word::null),
281  volumePhaseName_(word::null)
282 {
283  readCoeffs();
284  volumeAlpha();
285 }
286 
287 
288 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
289 
291 {}
292 
293 
294 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
295 
297 {
298  return true;
299 }
300 
301 
303 {
304  return wordList();
305 }
306 
307 
309 
310 
312 
313 
315 (
318 );
319 
320 
322 {
323  if (fvModel::read(dict))
324  {
325  readCoeffs();
326  return true;
327  }
328  else
329  {
330  return false;
331  }
332 }
333 
334 
335 // ************************************************************************* //
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
virtual ~volumeFractionSource()
Destructor.
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
virtual bool read(const dictionary &dict)
Read dictionary.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:158
static tmp< surfaceInterpolationScheme< Type > > scheme(const surfaceScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
Calculate the matrix for the laplacian of the field.
word member() const
Return member (name without the extension)
Definition: IOobject.C:336
Finite volume model abstract base class.
Definition: fvModel.H:55
word group() const
Return group (extension part of name)
Definition: IOobject.C:330
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
CompressibleMomentumTransportModel< dynamicTransportModel > momentumTransportModel
const dimensionSet dimTime
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static word groupName(Name name, const word &group)
FOR_ALL_FIELD_TYPES(IMPLEMENT_FV_MODEL_ADD_SUP, fv::volumeFractionSource)
This fvModel adds transport terms into the equations to account for the presence of a constant volume...
static const word null
An empty word.
Definition: word.H:77
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
Definition: fvModelM.H:33
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
const dimensionSet dimMass
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
volumeFractionSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
Calculate the matrix for the divergence of the given field and flux.
List< word > wordList
A List of words.
Definition: fileName.H:54
#define IMPLEMENT_FV_MODEL_ADD_RHO_SUP(Type, modelType)
Definition: fvModelM.H:51
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
const dimensionSet dimVolume
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.