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-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 "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  option,
44  volumeFractionSource,
45  dictionary
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 const Foam::volScalarField& Foam::fv::volumeFractionSource::alpha() const
54 {
55  const word alphaName = IOobject::groupName("alpha", phaseName_);
56 
57  if (!mesh_.foundObject<volScalarField>(alphaName))
58  {
59  volScalarField* alphaPtr =
60  new volScalarField
61  (
62  IOobject
63  (
64  alphaName,
65  mesh_.time().constant(),
66  mesh_,
69  ),
70  mesh_
71  );
72 
73  alphaPtr->store();
74  }
75 
76  return mesh_.lookupObject<volScalarField>(alphaName);
77 }
78 
79 
80 Foam::tmp<Foam::volScalarField> Foam::fv::volumeFractionSource::D
81 (
82  const label fieldi
83 ) const
84 {
85  const surfaceScalarField& phi =
86  mesh().lookupObject<surfaceScalarField>(phiName_);
87 
88  if (phi.dimensions() == dimVolume/dimTime)
89  {
90  const momentumTransportModel& turbulence =
91  mesh().lookupObject<momentumTransportModel>
92  (
93  momentumTransportModel::typeName
94  );
95 
96  return turbulence.nuEff();
97  }
98  else if (phi.dimensions() == dimMass/dimTime)
99  {
100  const thermophysicalTransportModel& ttm =
101  mesh().lookupObject<thermophysicalTransportModel>
102  (
103  thermophysicalTransportModel::typeName
104  );
105 
106  return
107  fieldNames_[fieldi] == ttm.thermo().T().name()
108  ? ttm.kappaEff()
109  : fieldNames_[fieldi] == ttm.thermo().he().name()
110  ? ttm.alphaEff()
111  : ttm.momentumTransport().muEff();
112  }
113  else
114  {
116  << "Dimensions of " << phi.name() << " not recognised"
117  << exit(FatalError);
118  return tmp<volScalarField>(nullptr);
119  }
120 }
121 
122 
123 template <class Type>
124 void Foam::fv::volumeFractionSource::addDivSup
125 (
126  fvMatrix<Type>& eqn,
127  const label fieldi
128 )
129 {
130  const surfaceScalarField& phi =
131  mesh().lookupObject<surfaceScalarField>(phiName_);
132 
133  const volScalarField AByB(this->alpha()/(1 - this->alpha()));
134 
135  eqn -= AByB*fvm::div(phi, eqn.psi());
136 }
137 
138 
139 void Foam::fv::volumeFractionSource::addUDivSup
140 (
141  fvMatrix<vector>& eqn,
142  const label fieldi
143 )
144 {
145  const surfaceScalarField& phi =
146  mesh().lookupObject<surfaceScalarField>(phiName_);
147 
148  const volScalarField AByB(this->alpha()/(1 - this->alpha()));
149 
150  const word scheme("div(" + phiName_ + "," + eqn.psi().name() + ")");
151 
152  eqn -= fvm::div(fvc::interpolate(AByB)*phi, eqn.psi(), scheme);
153 }
154 
155 
156 void Foam::fv::volumeFractionSource::addRhoDivSup
157 (
158  fvMatrix<scalar>& eqn,
159  const label fieldi
160 )
161 {
162  const surfaceScalarField& phi =
163  mesh().lookupObject<surfaceScalarField>(phiName_);
164 
165  const volScalarField AByB(this->alpha()/(1 - this->alpha()));
166 
167  eqn -= AByB*fvc::div(phi);
168 }
169 
170 
171 template <class Type, class AlphaFieldType>
172 void Foam::fv::volumeFractionSource::addLaplacianSup
173 (
174  const AlphaFieldType& alpha,
175  fvMatrix<Type>& eqn,
176  const label fieldi
177 )
178 {
179  const volScalarField B(1 - this->alpha());
180 
181  const volScalarField D(this->D(fieldi));
182 
183  const word scheme("laplacian(" + D.name() + "," + eqn.psi().name() + ")");
184 
185  eqn +=
186  fvm::laplacian(D, eqn.psi())
187  - 1/B*fvm::laplacian(B*D, eqn.psi(), scheme);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
194 (
195  const word& name,
196  const word& modelType,
197  const dictionary& dict,
198  const fvMesh& mesh
199 )
200 :
201  option(name, modelType, dict, mesh),
202  phaseName_(dict.lookup<word>("phase")),
203  phiName_("phi"),
204  rhoName_("rho"),
205  UName_("U")
206 {
207  read(dict);
208  alpha();
209 }
210 
211 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
212 
214 {}
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
220 (
221  fvMatrix<scalar>& eqn,
222  const label fieldi
223 )
224 {
225  if (fieldNames_[fieldi] == rhoName_)
226  {
227  addRhoDivSup(eqn, fieldi);
228  }
229  else
230  {
231  addDivSup(eqn, fieldi);
232  addLaplacianSup(geometricOneField(), eqn, fieldi);
233  }
234 }
235 
236 
238 (
239  fvMatrix<vector>& eqn,
240  const label fieldi
241 )
242 {
243  if (fieldNames_[fieldi] == UName_)
244  {
245  addUDivSup(eqn, fieldi);
246  }
247  else
248  {
249  addDivSup(eqn, fieldi);
250  addLaplacianSup(geometricOneField(), eqn, fieldi);
251  }
252 }
253 
254 
256 (
258  const label fieldi
259 )
260 {
261  addDivSup(eqn, fieldi);
262  addLaplacianSup(geometricOneField(), eqn, fieldi);
263 }
264 
265 
267 (
269  const label fieldi
270 )
271 {
272  addDivSup(eqn, fieldi);
273  addLaplacianSup(geometricOneField(), eqn, fieldi);
274 }
275 
276 
278 (
279  fvMatrix<tensor>& eqn,
280  const label fieldi
281 )
282 {
283  addDivSup(eqn, fieldi);
284  addLaplacianSup(geometricOneField(), eqn, fieldi);
285 }
286 
287 
289 (
290  const volScalarField& rho,
291  fvMatrix<scalar>& eqn,
292  const label fieldi
293 )
294 {
295  if (fieldNames_[fieldi] == rhoName_)
296  {
297  addRhoDivSup(eqn, fieldi);
298  }
299  else
300  {
301  addDivSup(eqn, fieldi);
302  addLaplacianSup(geometricOneField(), eqn, fieldi);
303  }
304 }
305 
306 
308 (
309  const volScalarField& rho,
310  fvMatrix<vector>& eqn,
311  const label fieldi
312 )
313 {
314  if (fieldNames_[fieldi] == UName_)
315  {
316  addUDivSup(eqn, fieldi);
317  }
318  else
319  {
320  addDivSup(eqn, fieldi);
321  addLaplacianSup(geometricOneField(), eqn, fieldi);
322  }
323 }
324 
325 
327 (
328  const volScalarField& rho,
330  const label fieldi
331 )
332 {
333  addDivSup(eqn, fieldi);
334  addLaplacianSup(geometricOneField(), eqn, fieldi);
335 }
336 
337 
339 (
340  const volScalarField& rho,
342  const label fieldi
343 )
344 {
345  addDivSup(eqn, fieldi);
346  addLaplacianSup(geometricOneField(), eqn, fieldi);
347 }
348 
349 
351 (
352  const volScalarField& rho,
353  fvMatrix<tensor>& eqn,
354  const label fieldi
355 )
356 {
357  addDivSup(eqn, fieldi);
358  addLaplacianSup(geometricOneField(), eqn, fieldi);
359 }
360 
361 
363 (
364  const volScalarField& alpha,
365  const volScalarField& rho,
366  fvMatrix<scalar>& eqn,
367  const label fieldi
368 )
369 {
370  if (fieldNames_[fieldi] == rhoName_)
371  {
372  addRhoDivSup(eqn, fieldi);
373  }
374  else
375  {
376  addDivSup(eqn, fieldi);
377  addLaplacianSup(alpha, eqn, fieldi);
378  }
379 }
380 
381 
383 (
384  const volScalarField& alpha,
385  const volScalarField& rho,
386  fvMatrix<vector>& eqn,
387  const label fieldi
388 )
389 {
390  if (fieldNames_[fieldi] == UName_)
391  {
392  addUDivSup(eqn, fieldi);
393  }
394  else
395  {
396  addDivSup(eqn, fieldi);
397  addLaplacianSup(alpha, eqn, fieldi);
398  }
399 }
400 
401 
403 (
404  const volScalarField& alpha,
405  const volScalarField& rho,
407  const label fieldi
408 )
409 {
410  addDivSup(eqn, fieldi);
411  addLaplacianSup(alpha, eqn, fieldi);
412 }
413 
414 
416 (
417  const volScalarField& alpha,
418  const volScalarField& rho,
420  const label fieldi
421 )
422 {
423  addDivSup(eqn, fieldi);
424  addLaplacianSup(alpha, eqn, fieldi);
425 }
426 
427 
429 (
430  const volScalarField& alpha,
431  const volScalarField& rho,
432  fvMatrix<tensor>& eqn,
433  const label fieldi
434 )
435 {
436  addDivSup(eqn, fieldi);
437  addLaplacianSup(alpha, eqn, fieldi);
438 }
439 
440 
442 {
443  if (option::read(dict))
444  {
445  if (coeffs_.found("fields"))
446  {
447  coeffs_.lookup("fields") >> fieldNames_;
448  }
449 
450  applied_.setSize(fieldNames_.size(), false);
451 
452  dict.readIfPresent("phi", phiName_);
453 
454  dict.readIfPresent("rho", rhoName_);
455 
456  dict.readIfPresent("U", UName_);
457 
458  return true;
459  }
460  else
461  {
462  return false;
463  }
464 }
465 
466 
467 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
virtual ~volumeFractionSource()
Destructor.
virtual void addSup(fvMatrix< scalar > &, const label)
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:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool read(const dictionary &dict)
Read dictionary.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
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.
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
CompressibleMomentumTransportModel< fluidThermo > momentumTransportModel
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static word groupName(Name name, const word &group)
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 &dict)
Read source dictionary.
Definition: fvOptionIO.C:53
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
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
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.
Calculate the matrix for the divergence of the given field and flux.
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
Finite volume options abstract base class. Provides a base set of controls, e.g.: ...
Definition: fvOption.H:66