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-2023 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 {
42  (
43  fvModel,
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  {
103  mesh().lookupType<momentumTransportModel>();
104 
105  return turbulence.nuEff();
106  }
107  else if (phi.dimensions() == dimMass/dimTime)
108  {
109  const fluidThermophysicalTransportModel& ttm =
110  mesh().lookupType<fluidThermophysicalTransportModel>();
111 
112  return
113  fieldName == ttm.thermo().T().name()
114  ? ttm.kappaEff()
115  : fieldName == ttm.thermo().he().name()
116  ? ttm.kappaEff()/ttm.thermo().Cpv()
117  : ttm.momentumTransport().rho()*ttm.momentumTransport().nuEff();
118  }
119  else
120  {
122  << "Dimensions of " << phi.name() << " not recognised"
123  << exit(FatalError);
124 
125  return tmp<volScalarField>(nullptr);
126  }
127 }
128 
129 
130 template <class Type, class AlphaFieldType>
131 void Foam::fv::volumeFractionSource::addGeneralSup
132 (
133  const AlphaFieldType& alpha,
134  fvMatrix<Type>& eqn,
135  const word& fieldName
136 ) const
137 {
138  const word phiName =
139  IOobject::groupName(phiName_, IOobject::group(fieldName));
140  const surfaceScalarField& phi =
141  mesh().lookupObject<surfaceScalarField>(phiName);
142 
143  const volScalarField B(1 - volumeAlpha());
144  const volScalarField AByB(volumeAlpha()/B);
145  const volScalarField D(this->D(fieldName));
146 
147  // Divergence term
148  const word divScheme = "div(" + phiName + "," + eqn.psi().name() + ")";
149  eqn -= AByB*fvm::div(phi, eqn.psi(), divScheme);
150 
151  // Laplacian term
152  const word laplacianScheme =
153  "laplacian(" + D.name() + "," + eqn.psi().name() + ")";
154  eqn +=
155  fvm::laplacian(D, eqn.psi())
156  - 1/B*fvm::laplacian(B*D, eqn.psi(), laplacianScheme);
157 }
158 
159 
160 template<class Type, class AlphaFieldType>
161 void Foam::fv::volumeFractionSource::addAlphaSupType
162 (
163  const AlphaFieldType& alpha,
164  fvMatrix<Type>& eqn,
165  const word& fieldName
166 ) const
167 {
168  addGeneralSup(alpha, eqn, fieldName);
169 }
170 
171 
172 template<class AlphaFieldType>
173 void Foam::fv::volumeFractionSource::addAlphaSupType
174 (
175  const AlphaFieldType& alpha,
176  fvMatrix<scalar>& eqn,
177  const word& fieldName
178 ) const
179 {
180  if (IOobject::member(fieldName) == rhoName_)
181  {
182  const word phiName =
183  IOobject::groupName(phiName_, IOobject::group(fieldName));
184  const surfaceScalarField& phi =
185  mesh().lookupObject<surfaceScalarField>(phiName);
186 
187  const volScalarField AByB(volumeAlpha()/(1 - volumeAlpha()));
188 
189  eqn -= AByB*fvc::div(phi);
190  }
191  else
192  {
193  addGeneralSup(alpha, eqn, fieldName);
194  }
195 }
196 
197 
198 template<class AlphaFieldType>
199 void Foam::fv::volumeFractionSource::addAlphaSupType
200 (
201  const AlphaFieldType& alpha,
202  fvMatrix<vector>& eqn,
203  const word& fieldName
204 ) const
205 {
206  if (IOobject::member(fieldName) == UName_)
207  {
208  const word phiName =
209  IOobject::groupName(phiName_, IOobject::group(fieldName));
210  const surfaceScalarField& phi =
211  mesh().lookupObject<surfaceScalarField>(phiName);
212 
213  const volScalarField AByB(volumeAlpha()/(1 - volumeAlpha()));
214 
215  const word scheme("div(" + phiName + "," + eqn.psi().name() + ")");
216 
217  eqn -= fvm::div(fvc::interpolate(AByB)*phi, eqn.psi(), scheme);
218  }
219  else
220  {
221  addGeneralSup(alpha, eqn, fieldName);
222  }
223 }
224 
225 
226 template<class Type>
227 void Foam::fv::volumeFractionSource::addSupType
228 (
229  fvMatrix<Type>& eqn,
230  const word& fieldName
231 ) const
232 {
233  addAlphaSupType(geometricOneField(), eqn, fieldName);
234 }
235 
236 
237 template<class Type>
238 void Foam::fv::volumeFractionSource::addSupType
239 (
240  const volScalarField& rho,
241  fvMatrix<Type>& eqn,
242  const word& fieldName
243 ) const
244 {
245  addAlphaSupType(geometricOneField(), eqn, fieldName);
246 }
247 
248 
249 template<class Type>
250 void Foam::fv::volumeFractionSource::addSupType
251 (
252  const volScalarField& alpha,
253  const volScalarField& rho,
254  fvMatrix<Type>& eqn,
255  const word& fieldName
256 ) const
257 {
258  addAlphaSupType(alpha, eqn, fieldName);
259 }
260 
261 
262 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
263 
265 (
266  const word& name,
267  const word& modelType,
268  const fvMesh& mesh,
269  const dictionary& dict
270 )
271 :
272  fvModel(name, modelType, mesh, dict),
273  phiName_(word::null),
274  rhoName_(word::null),
275  UName_(word::null),
276  volumePhaseName_(word::null)
277 {
278  readCoeffs();
279  volumeAlpha();
280 }
281 
282 
283 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
284 
286 {}
287 
288 
289 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
290 
292 {
293  return true;
294 }
295 
296 
298 {
299  return wordList();
300 }
301 
302 
304 
305 
307 
308 
310 (
312  fv::volumeFractionSource
313 );
314 
315 
317 {
318  return true;
319 }
320 
321 
323 {}
324 
325 
327 {}
328 
329 
331 {}
332 
333 
335 {
336  if (fvModel::read(dict))
337  {
338  readCoeffs();
339  return true;
340  }
341  else
342  {
343  return false;
344  }
345 }
346 
347 
348 // ************************************************************************* //
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
word member() const
Return member (name without the extension)
Definition: IOobject.C:330
static word groupName(Name name, const word &group)
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
This fvModel adds transport terms into the equations to account for the presence of a constant volume...
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual ~volumeFractionSource()
Destructor.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
volumeFractionSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define IMPLEMENT_FV_MODEL_ADD_RHO_SUP(Type, modelType)
Definition: fvModelM.H:51
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
Definition: fvModelM.H:33
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
Calculate the divergence of the given field.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
compressibleMomentumTransportModel momentumTransportModel
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
static tmp< surfaceInterpolationScheme< Type > > scheme(const surfaceScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
SurfaceField< scalar > surfaceScalarField
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
const dimensionSet dimTime
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
const dimensionSet dimMass
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList fv(nPoints)
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))