volumeSource.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) 2021-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 "volumeSource.H"
27 #include "fvMatrices.H"
28 #include "basicThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::fv::volumeSource::readCoeffs()
46 {
47  alphaName_ =
49  ? word::null
50  : coeffs().lookupOrDefault<word>
51  (
52  "alpha",
53  IOobject::groupName("alpha", phaseName())
54  );
55 
56  setPtr_->read(coeffs());
57 
58  volumetricFlowRate_.reset
59  (
61  (
62  "volumetricFlowRate",
63  mesh().time().userUnits(),
65  coeffs()
66  ).ptr()
67  );
68 }
69 
70 
71 template<class Type>
72 void Foam::fv::volumeSource::addSupType
73 (
74  const VolField<Type>& field,
75  fvMatrix<Type>& eqn
76 ) const
77 {
79  << "field=" << field.name()
80  << ", eqnField=" << eqn.psi().name() << endl;
81 
82  // Single-phase property equation
83  if (phaseName() == word::null && field.group() == word::null)
84  {
85  fvTotalSource::addSupType(field, eqn);
86  }
87  // Multiphase volume-weighted mixture property equation (e.g., a turbulence
88  // equation if running standard incompressible transport modelling in the
89  // incompressibleVoF solver)
90  else if (phaseName() != word::null && field.group() == word::null)
91  {
92  fvTotalSource::addSupType(field, eqn);
93  }
94  // Not recognised. Fail.
95  else
96  {
97  const volScalarField& null = NullObjectRef<volScalarField>();
98  addSupType(null, null, field, eqn);
99  }
100 }
101 
102 
103 void Foam::fv::volumeSource::addSupType
104 (
105  const volScalarField& alphaOrField,
106  fvMatrix<scalar>& eqn
107 ) const
108 {
110  << "alphaOrField=" << alphaOrField.name()
111  << ", eqnField=" << eqn.psi().name() << endl;
112 
113  // Multiphase continuity equation
114  if (phaseName() != word::null && alphaOrField.name() == alphaName_)
115  {
117  }
118  // Try the general type method
119  else
120  {
121  addSupType<scalar>(alphaOrField, eqn);
122  }
123 }
124 
125 
126 template<class Type>
127 void Foam::fv::volumeSource::addSupType
128 (
129  const volScalarField& alphaOrRho,
130  const VolField<Type>& field,
131  fvMatrix<Type>& eqn
132 ) const
133 {
135  << "alphaOrRho=" << alphaOrRho.name()
136  << ", field=" << field.name()
137  << ", eqnField=" << eqn.psi().name() << endl;
138 
139  // Multiphase property equation (e.g., a turbulence equation if running
140  // two-phase transport modelling in the incompressibleVoF solver)
141  if (phaseName() != word::null && alphaOrRho.name() == alphaName_)
142  {
143  fvTotalSource::addSupType(field, eqn);
144  }
145  // Multiphase mass-weighted mixture property equation (e.g., the momentum
146  // equation in the incompressibleVoF solver)
147  else if
148  (
149  phaseName() != word::null
150  && alphaOrRho.group() == word::null
151  && alphaOrRho.dimensions() == dimDensity
152  && field.group() == word::null
153  )
154  {
155  // First we construct the volumetric source...
156  fvMatrix<Type> volEqn(eqn.psi(), eqn.dimensions()/dimDensity);
157  fvTotalSource::addSupType(field, volEqn);
158 
159  // Then, to apply it to the mixture equation, we need to multiply by
160  // the density of the phase of which this is a source. There is no
161  // solver-agnostic interface, at present, that lets us obtain this
162  // density. So, we read it from the physical properties file. This is
163  // clunky, but it should work in all circumstances. This is what the
164  // clouds fvModel does,
165  const dimensionedScalar rhoi
166  (
167  "rho",
168  dimDensity,
169  mesh().lookupObject<IOdictionary>
170  (
172  (
173  physicalProperties::typeName,
174  phaseName()
175  )
176  )
177  );
178  eqn += rhoi*volEqn;
179  }
180  // Not recognised. Fail.
181  else
182  {
183  const volScalarField& null = NullObjectRef<volScalarField>();
184  addSupType(null, alphaOrRho, field, eqn);
185  }
186 }
187 
188 
189 template<class Type>
190 void Foam::fv::volumeSource::addSupType
191 (
192  const volScalarField& alpha,
193  const volScalarField& rho,
194  const VolField<Type>& field,
195  fvMatrix<Type>& eqn
196 ) const
197 {
199  << "alpha=" << (isNull(alpha) ? word::null : alpha.name())
200  << ", rho=" << (isNull(rho) ? word::null : rho.name())
201  << ", field=" << field.name()
202  << ", eqnField=" << eqn.psi().name() << endl;
203 
205  << "Cannot add a volume source for field " << field.name()
206  << " to equation for " << eqn.psi().name() << " because this field's "
207  << "equation was not recognised as being in volume-conservative form"
208  << exit(FatalError);
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
213 
215 (
216  const word& name,
217  const word& modelType,
218  const fvMesh& mesh,
219  const dictionary& dict
220 )
221 :
222  fvTotalSource(name, modelType, mesh, dict),
223  alphaName_(),
224  setPtr_(new fvCellSet(mesh)),
225  volumetricFlowRate_()
226 {
227  readCoeffs();
228 }
229 
230 
231 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
232 
234 {
235  return setPtr_->cells();
236 }
237 
238 
240 {
241  return setPtr_->nCells();
242 }
243 
244 
245 Foam::scalar Foam::fv::volumeSource::V() const
246 {
247  return setPtr_->V();
248 }
249 
250 
252 {
253  return
255  (
257  volumetricFlowRate_->value(mesh().time().value())
258  );
259 }
260 
261 
263 {
265  << "eqnField=" << eqn.psi().name() << endl;
266 
267  // Single-phase continuity equation
269 }
270 
271 
273 
274 
276 
277 
279 (
281  fv::volumeSource
282 )
283 
284 
285 bool Foam::fv::volumeSource::movePoints()
286 {
287  setPtr_->movePoints();
288  return true;
289 }
290 
291 
293 {
294  setPtr_->topoChange(map);
295 }
296 
297 
299 {
300  setPtr_->mapMesh(map);
301 }
302 
303 
305 {
306  setPtr_->distribute(map);
307 }
308 
309 
311 {
313  {
314  readCoeffs();
315  return true;
316  }
317  else
318  {
319  return false;
320  }
321 }
322 
323 
324 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< scalar > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
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.
General run-time selected cell set selection class for fvMesh.
Definition: fvCellSet.H:88
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
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
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
Base class for sources which are specified as a total value (e.g., volume or mass flow rate),...
Definition: fvTotalSource.H:55
void addSupType(const VolField< Type > &field, fvMatrix< Type > &eqn) const
Add a source term to an equation.
const word & phaseName() const
Return the phase name.
virtual bool read(const dictionary &dict)
Read source dictionary.
void addSource(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: fvTotalSource.C:48
This fvModel applies a volume source to the continuity equation and to all field equations....
Definition: volumeSource.H:90
virtual scalar V() const
Return the volume of cells that the source applies to.
Definition: volumeSource.C:245
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: volumeSource.C:292
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: volumeSource.C:304
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: volumeSource.C:310
virtual labelUList cells() const
Return the cells that the source applies to.
Definition: volumeSource.C:233
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: volumeSource.C:298
volumeSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: volumeSource.C:215
virtual label nCells() const
Return the number of cells that the source applies to.
Definition: volumeSource.C:239
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: volumeSource.C:262
virtual dimensionedScalar S() const
Return the source value.
Definition: volumeSource.C:251
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 handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
A special matrix type and solver, designed for finite volume solutions of scalar equations.
#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
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define DebugInFunction
Report an information message using Foam::Info.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict