massTransfer.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 "massTransfer.H"
27 #include "fvMatrices.H"
28 #include "physicalProperties.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::fv::massTransfer::readCoeffs()
45 {
46  if
47  (
48  phaseNames_ != lookupPhaseNames()
49  || alphaNames_ != lookupPhaseFieldNames("alpha")
50  || rhoNames_ != lookupPhaseFieldNames("rho")
51  )
52  {
54  << "Cannot change the phases of a " << typeName << " model "
55  << "at run time" << exit(FatalIOError);
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
61 
63 {
64  return coeffs().lookup<Pair<word>>("phases");
65 }
66 
67 
70 {
71  return
72  coeffs().lookupOrDefault<Pair<word>>
73  (
74  name + "s",
76  (
77  IOobject::groupName(name, phaseNames_.first()),
78  IOobject::groupName(name, phaseNames_.second())
79  )
80  );
81 }
82 
83 
85 (
86  const label i
87 ) const
88 {
89  // Compressible case. Lookup the field.
90  if (mesh().foundObject<volScalarField::Internal>(rhoNames()[i]))
91  {
92  return mesh().lookupObject<volScalarField::Internal>(rhoNames()[i]);
93  }
94 
95  // Incompressible case. Read from the physical properties dictionary.
96  const word physicalPropertiesName =
97  IOobject::groupName(physicalProperties::typeName, phaseNames()[i]);
98 
99  if (mesh().foundObject<IOdictionary>(physicalPropertiesName))
100  {
102  mesh().lookupObject<IOdictionary>(physicalPropertiesName);
103 
104  if (physicalProperties.found("rho"))
105  {
106  return
108  (
109  rhoNames()[i],
110  mesh(),
112  );
113  }
114  }
115 
116  // Fail.
118  << "Could not determine the density " << rhoNames()[i]
119  << " for phase " << phaseNames()[i]
120  << exit(FatalError);
121  return tmp<volScalarField::Internal>(nullptr);
122 }
123 
124 
126 (
127  const volScalarField& alphaOrField,
128  fvMatrix<scalar>& eqn
129 ) const
130 {
132  << "alphaOrField=" << alphaOrField.name()
133  << ", eqnField=" << eqn.psi().name() << endl;
134 
135  const label i = index(alphaNames(), alphaOrField.name());
136 
137  // Incompressible continuity equation
138  if (i != -1)
139  {
140  eqn += S(alphaOrField.name())/rho(i);
141  }
142  // Not recognised. Fail.
143  else
144  {
145  addSupType<scalar>(alphaOrField, eqn);
146  }
147 }
148 
149 
151 (
152  const volScalarField& alphaOrRho,
153  const volScalarField& rhoOrField,
154  fvMatrix<scalar>& eqn
155 ) const
156 {
158  << "alphaOrRho=" << alphaOrRho.name()
159  << ", rhoOrField=" << rhoOrField.name()
160  << ", eqnField=" << eqn.psi().name() << endl;
161 
162  // Compressible continuity equation
163  if
164  (
165  index(alphaNames(), alphaOrRho.name()) != -1
166  && index(rhoNames(), rhoOrField.name()) != -1
167  )
168  {
169  eqn += S(alphaOrRho.name());
170  }
171  // Mixture property equation
172  else
173  {
174  addSupType<scalar>(alphaOrRho, rhoOrField, eqn);
175  }
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
180 
182 (
183  const word& name,
184  const word& modelType,
185  const fvMesh& mesh,
186  const dictionary& dict
187 )
188 :
189  fvSpecificSource(name, modelType, mesh, dict),
190  phaseNames_(lookupPhaseNames()),
191  alphaNames_(lookupPhaseFieldNames("alpha")),
192  rhoNames_(lookupPhaseFieldNames("rho"))
193 {
194  readCoeffs();
195 }
196 
197 
198 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
199 
200 bool Foam::fv::massTransfer::addsSupToField(const word& fieldName) const
201 {
202  const word group = IOobject::group(fieldName);
203 
204  return
205  group == word::null
206  || group == phaseNames_.first()
207  || group == phaseNames_.second();
208 }
209 
210 
212 Foam::fv::massTransfer::S(const word& fieldName) const
213 {
214  return sign(phaseNames(), IOobject::group(fieldName))*mDot();
215 }
216 
217 
219 {
221  << "eqnField=" << eqn.psi().name() << endl;
222 
224  << "Field-less mass transfers are not possible"
225  << exit(FatalError);
226 }
227 
228 
230 
231 
233 
234 
236 (
238  fv::massTransfer
239 );
240 
241 
243 {
245  {
246  readCoeffs();
247  return true;
248  }
249  else
250  {
251  return false;
252  }
253 }
254 
255 
256 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
const word & name() const
Return name.
Definition: IOobject.H:310
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
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
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
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
Base class for sources which are specified as a specific value (e.g., mass flow rate per unit volume)...
Base class for mass transfers between phases.
Definition: massTransfer.H:55
void addSupType(const VolField< Type > &field, fvMatrix< Type > &eqn) const
Add a source term to an equation.
const Pair< word > lookupPhaseFieldNames(const word &name) const
Lookup the phase field names.
Definition: massTransfer.C:69
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massTransfer.C:242
tmp< volScalarField::Internal > rho(const label i) const
Return the density.
Definition: massTransfer.C:85
massTransfer(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: massTransfer.C:182
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:218
const Pair< word > lookupPhaseNames() const
Lookup the phase names.
Definition: massTransfer.C:62
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: massTransfer.C:200
virtual tmp< DimensionedField< scalar, volMesh > > S(const word &fieldName) const
Return the source value.
Definition: massTransfer.C:212
A base class for physical properties.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#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
#define DebugInFunction
Report an information message using Foam::Info.
const char *const group
Group name for atomic constants.
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar sign(const dimensionedScalar &ds)
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 dimDensity
IOerror FatalIOError
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict