massSourceBase.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 "massSourceBase.H"
27 #include "fvMatrices.H"
29 
30 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::fv::massSourceBase::readCoeffs()
44 {
45  rhoName_ =
46  coeffs().lookupOrDefault<word>
47  (
48  "rho",
50  );
51 }
52 
53 
54 template<class Type>
55 void Foam::fv::massSourceBase::addSupType
56 (
57  const VolField<Type>& field,
58  fvMatrix<Type>& eqn
59 ) const
60 {
62  << "field=" << field.name()
63  << ", eqnField=" << eqn.psi().name() << endl;
64 
66  << "Cannot add a mass source for field " << field.name()
67  << " to equation for " << eqn.psi().name() << " because this field's "
68  << "equation was not recognised as being in mass-conservative form"
69  << exit(FatalError);
70 }
71 
72 
73 void Foam::fv::massSourceBase::addSupType
74 (
75  const volScalarField& rhoOrField,
76  fvMatrix<scalar>& eqn
77 ) const
78 {
80  << "rhoOrField=" << rhoOrField.name()
81  << ", eqnField=" << eqn.psi().name() << endl;
82 
83  // Single-phase continuity equation
84  if (phaseName() == word::null && rhoOrField.name() == rhoName_)
85  {
87  }
88  // Not recognised. Fail.
89  else
90  {
91  addSupType<scalar>(rhoOrField, eqn);
92  }
93 }
94 
95 
96 template<class Type>
97 void Foam::fv::massSourceBase::addSupType
98 (
99  const volScalarField& rho,
100  const VolField<Type>& field,
101  fvMatrix<Type>& eqn
102 ) const
103 {
105  << "rho=" << rho.name()
106  << ", field=" << field.name()
107  << ", eqnField=" << eqn.psi().name() << endl;
108 
109  // Single-phase property equation
110  if (phaseName() == word::null && rho.name() == rhoName_)
111  {
112  fvTotalSource::addSupType(rho, field, eqn);
113  }
114  // Multiphase mass-weighted mixture property equation
115  else if
116  (
117  phaseName() != word::null
118  && rho.group() == word::null
119  && rho.dimensions() == dimDensity
120  && field.group() == word::null
121  )
122  {
123  fvTotalSource::addSupType(rho, field, eqn);
124  }
125  // Not recognised. Fail.
126  else
127  {
128  addSupType<Type>(field, eqn);
129  }
130 }
131 
132 
133 void Foam::fv::massSourceBase::addSupType
134 (
135  const volScalarField& alphaOrRho,
136  const volScalarField& rhoOrField,
137  fvMatrix<scalar>& eqn
138 ) const
139 {
141  << "alphaOrRho=" << alphaOrRho.name()
142  << ", rhoOrField=" << rhoOrField.name()
143  << ", eqnField=" << eqn.psi().name() << endl;
144 
145  // Multiphase continuity equation
146  if (rhoOrField.name() == rhoName_)
147  {
149  }
150  // Try the general type method
151  else
152  {
153  addSupType<scalar>(alphaOrRho, rhoOrField, eqn);
154  }
155 }
156 
157 
158 template<class Type>
159 void Foam::fv::massSourceBase::addSupType
160 (
161  const volScalarField& alpha,
162  const volScalarField& rho,
163  const VolField<Type>& field,
164  fvMatrix<Type>& eqn
165 ) const
166 {
168  << "alpha=" << alpha.name()
169  << ", rho=" << rho.name()
170  << ", field=" << field.name()
171  << ", eqnField=" << eqn.psi().name() << endl;
172 
173  // Multiphase property equation
174  fvTotalSource::addSupType(alpha, rho, field, eqn);
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
179 
181 (
182  const word& name,
183  const word& modelType,
184  const fvMesh& mesh,
185  const dictionary& dict
186 )
187 :
188  fvTotalSource(name, modelType, mesh, dict),
189  rhoName_()
190 {
191  readCoeffs();
192 }
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
198 {
200  << "eqnField=" << eqn.psi().name() << endl;
201 
203  << "Field-less mass sources are not possible"
204  << exit(FatalError);
205 }
206 
207 
209 
210 
212 
213 
215 (
217  fv::massSourceBase
218 );
219 
220 
222 {
224  {
225  readCoeffs();
226  return true;
227  }
228  else
229  {
230  return false;
231  }
232 }
233 
234 
235 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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.
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
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
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
Base class for mass source models.
virtual bool read(const dictionary &dict)
Read source dictionary.
massSourceBase(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
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.
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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 dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
labelList fv(nPoints)
dictionary dict