volumeBlockage.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-2024 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 "volumeBlockage.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 {
43  (
44  fvModel,
46  dictionary,
47  volumeFractionSource,
48  "volumeFractionSource"
49  );
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void Foam::fv::volumeBlockage::readCoeffs()
57 {
58  phiName_ = coeffs().lookupOrDefault<word>("phi", "phi");
59  rhoName_ = coeffs().lookupOrDefault<word>("rho", "rho");
60  UName_ = coeffs().lookupOrDefault<word>("U", "U");
61 
62  volumePhaseName_ = coeffs().lookup<word>("volumePhase");
63 }
64 
65 
66 const Foam::volScalarField& Foam::fv::volumeBlockage::volumeAlpha() const
67 {
68  const word alphaName = IOobject::groupName("alpha", volumePhaseName_);
69 
70  if (!mesh().foundObject<volScalarField>(alphaName))
71  {
72  volScalarField* alphaPtr =
73  new volScalarField
74  (
75  IOobject
76  (
77  alphaName,
78  mesh().time().constant(),
79  mesh(),
82  ),
83  mesh()
84  );
85 
86  alphaPtr->store();
87  }
88 
89  return mesh().lookupObject<volScalarField>(alphaName);
90 }
91 
92 
93 Foam::tmp<Foam::volScalarField> Foam::fv::volumeBlockage::D
94 (
95  const word& fieldName
96 ) const
97 {
98  const word phiName =
99  IOobject::groupName(phiName_, IOobject::group(fieldName));
100  const surfaceScalarField& phi =
101  mesh().lookupObject<surfaceScalarField>(phiName);
102 
103  if (phi.dimensions() == dimVolume/dimTime)
104  {
106  mesh().lookupType<momentumTransportModel>();
107 
108  return turbulence.nuEff();
109  }
110  else if (phi.dimensions() == dimMass/dimTime)
111  {
112  const fluidThermophysicalTransportModel& ttm =
113  mesh().lookupType<fluidThermophysicalTransportModel>();
114 
115  return
116  fieldName == ttm.thermo().T().name()
117  ? ttm.kappaEff()
118  : fieldName == ttm.thermo().he().name()
119  ? ttm.kappaEff()/ttm.thermo().Cpv()
120  : ttm.momentumTransport().rho()*ttm.momentumTransport().nuEff();
121  }
122  else
123  {
125  << "Dimensions of " << phi.name() << " not recognised"
126  << exit(FatalError);
127 
128  return tmp<volScalarField>(nullptr);
129  }
130 }
131 
132 
133 template <class Type, class AlphaFieldType>
134 void Foam::fv::volumeBlockage::addGeneralSupType
135 (
136  const AlphaFieldType& alpha,
137  fvMatrix<Type>& eqn
138 ) const
139 {
140  const word phiName =
141  IOobject::groupName(phiName_, IOobject::group(eqn.psi().name()));
142  const surfaceScalarField& phi =
143  mesh().lookupObject<surfaceScalarField>(phiName);
144 
145  const volScalarField B(1 - volumeAlpha());
146  const volScalarField AByB(volumeAlpha()/B);
147  const volScalarField D(this->D(eqn.psi().name()));
148 
149  // Divergence term
150  const word divScheme = "div(" + phiName + "," + eqn.psi().name() + ")";
151  eqn -= AByB*fvm::div(phi, eqn.psi(), divScheme);
152 
153  // Laplacian term
154  const word laplacianScheme =
155  "laplacian(" + D.name() + "," + eqn.psi().name() + ")";
156  eqn +=
157  fvm::laplacian(D, eqn.psi())
158  - 1/B*fvm::laplacian(B*D, eqn.psi(), laplacianScheme);
159 }
160 
161 
162 template<class Type, class AlphaFieldType>
163 void Foam::fv::volumeBlockage::addAlphaSupType
164 (
165  const AlphaFieldType& alpha,
166  const VolField<Type>& field,
167  fvMatrix<Type>& eqn
168 ) const
169 {
170  addGeneralSupType(alpha, eqn);
171 }
172 
173 
174 template<class AlphaFieldType>
175 void Foam::fv::volumeBlockage::addAlphaSupType
176 (
177  const AlphaFieldType& alpha,
178  const volScalarField& field,
179  fvMatrix<scalar>& eqn
180 ) const
181 {
182  if (IOobject::member(field.name()) == rhoName_)
183  {
184  const word phiName =
185  IOobject::groupName(phiName_, IOobject::group(field.name()));
186  const surfaceScalarField& phi =
187  mesh().lookupObject<surfaceScalarField>(phiName);
188 
189  const volScalarField AByB(volumeAlpha()/(1 - volumeAlpha()));
190 
191  eqn -= AByB*fvc::div(phi);
192  }
193  else
194  {
195  addGeneralSupType(alpha, eqn);
196  }
197 }
198 
199 
200 template<class AlphaFieldType>
201 void Foam::fv::volumeBlockage::addAlphaSupType
202 (
203  const AlphaFieldType& alpha,
204  const volVectorField& field,
205  fvMatrix<vector>& eqn
206 ) const
207 {
208  if (IOobject::member(field.name()) == UName_)
209  {
210  const word phiName =
211  IOobject::groupName(phiName_, IOobject::group(field.name()));
212  const surfaceScalarField& phi =
213  mesh().lookupObject<surfaceScalarField>(phiName);
214 
215  const volScalarField AByB(volumeAlpha()/(1 - volumeAlpha()));
216 
217  const word scheme("div(" + phiName + "," + eqn.psi().name() + ")");
218 
219  eqn -= fvm::div(fvc::interpolate(AByB)*phi, eqn.psi(), scheme);
220  }
221  else
222  {
223  addGeneralSupType(alpha, eqn);
224  }
225 }
226 
227 
228 template<class Type>
229 void Foam::fv::volumeBlockage::addSupType
230 (
231  const VolField<Type>& field,
232  fvMatrix<Type>& eqn
233 ) const
234 {
235  addAlphaSupType(geometricOneField(), field, eqn);
236 }
237 
238 
239 template<class Type>
240 void Foam::fv::volumeBlockage::addSupType
241 (
242  const volScalarField& rho,
243  const VolField<Type>& field,
244  fvMatrix<Type>& eqn
245 ) const
246 {
247  addAlphaSupType(geometricOneField(), field, eqn);
248 }
249 
250 
251 template<class Type>
252 void Foam::fv::volumeBlockage::addSupType
253 (
254  const volScalarField& alpha,
255  const volScalarField& rho,
256  const VolField<Type>& field,
257  fvMatrix<Type>& eqn
258 ) const
259 {
260  addAlphaSupType(alpha, field, eqn);
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265 
267 (
268  const word& name,
269  const word& modelType,
270  const fvMesh& mesh,
271  const dictionary& dict
272 )
273 :
274  fvModel(name, modelType, mesh, dict),
275  phiName_(word::null),
276  rhoName_(word::null),
277  UName_(word::null),
278  volumePhaseName_(word::null)
279 {
280  readCoeffs();
281  volumeAlpha();
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
286 
288 {}
289 
290 
291 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
292 
293 bool Foam::fv::volumeBlockage::addsSupToField(const word& fieldName) const
294 {
295  return true;
296 }
297 
298 
300 {
301  return wordList();
302 }
303 
304 
306 
307 
309 
310 
312 (
314  fv::volumeBlockage
315 )
316 
317 
318 bool Foam::fv::volumeBlockage::movePoints()
319 {
320  return true;
321 }
322 
323 
325 {}
326 
327 
329 {}
330 
331 
333 {}
334 
335 
337 {
338  if (fvModel::read(dict))
339  {
340  readCoeffs();
341  return true;
342  }
343  else
344  {
345  return false;
346  }
347 }
348 
349 
350 // ************************************************************************* //
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:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
This fvModel adds transport terms into the equations to account for the presence of a constant volume...
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
volumeBlockage(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
virtual ~volumeBlockage()
Destructor.
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:334
#define IMPLEMENT_FV_MODEL_ADD_FIELD_SUP(Type, modelType)
Definition: fvModelM.H:33
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_FIELD_SUP(Type, modelType)
Definition: fvModelM.H:71
#define IMPLEMENT_FV_MODEL_ADD_RHO_FIELD_SUP(Type, modelType)
Definition: fvModelM.H:51
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)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
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
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimTime
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
const dimensionSet dimMass
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
labelList fv(nPoints)
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))