fvTotalSource.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-2025 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 "fvTotalSource.H"
27 #include "fvMatrix.H"
28 
29 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34 }
35 
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 void Foam::fvTotalSource::readCoeffs(const dictionary& dict)
40 {
41  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
42 }
43 
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
50  << "eqnField=" << eqn.psi().name() << endl;
51 
52  const labelList& cells = this->zone();
53  const scalar V = this->V();
54  const dimensionedScalar S = this->S();
55 
56  // Check the dimensions
57  eqn.dimensions() = S.dimensions();
58 
59  // Apply the source
60  scalarField& eqnSource = eqn.source();
61  forAll(cells, i)
62  {
63  const scalar f = mesh().V()[cells[i]]/V;
64  eqnSource[cells[i]] -= f*S.value();
65  }
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const word& name,
74  const word& modelType,
75  const fvMesh& mesh,
76  const dictionary& dict
77 )
78 :
79  fvSource(name, modelType, mesh, dict),
80  phaseName_()
81 {
82  readCoeffs(coeffs(dict));
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
94 bool Foam::fvTotalSource::addsSupToField(const word& fieldName) const
95 {
96  const word group = IOobject::group(fieldName);
97 
98  return group == word::null || group == phaseName_;
99 }
100 
101 
103 {
104  if (fvModel::read(dict))
105  {
106  readCoeffs(coeffs(dict));
107  return true;
108  }
109  else
110  {
111  return false;
112  }
113 }
114 
115 
116 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
word group() const
Return group (extension part of name)
Definition: IOobject.C:321
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
VolField< Type > & psi()
Definition: fvMatrix.H:289
const dimensionSet & dimensions() const
Definition: fvMatrix.H:302
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
Base class for finite volume sources.
Definition: fvSource.H:52
Base class for sources which are specified as a total value (e.g., volume or mass flow rate),...
Definition: fvTotalSource.H:54
fvTotalSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: fvTotalSource.C:72
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual ~fvTotalSource()
Destructor.
Definition: fvTotalSource.C:88
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: fvTotalSource.C:94
void addSource(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: fvTotalSource.C:47
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const cellShapeList & cells
#define DebugInFunction
Report an information message using Foam::Info.
const char *const group
Group name for atomic constants.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
labelList f(nPoints)
dictionary dict