fvTotalSource.H
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-2023 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 Class
25  Foam::fvTotalSource
26 
27 Description
28  Base class for sources which are specified as a total value (e.g., volume
29  or mass flow rate), rather than a specific value (e.g., mass flow rate per
30  unit volume).
31 
32 SourceFiles
33  fvTotalSource.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef fvTotalSource_H
38 #define fvTotalSource_H
39 
40 #include "fvSource.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 class fvCellSet;
48 
49 /*---------------------------------------------------------------------------*\
50  Class fvTotalSource Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 class fvTotalSource
54 :
55  public fvSource
56 {
57 private:
58 
59  // Private Data
60 
61  //- Name of the phase
62  word phaseName_;
63 
64 
65  // Private Member Functions
66 
67  //- Non-virtual read
68  void readCoeffs();
69 
70 
71 protected:
72 
73  // Protected Member Functions
74 
75  // Sources
76 
77  //- Add a source term to a field-less proxy equation
78  void addSource(fvMatrix<scalar>& eqn) const;
79 
80  //- Add a source term to an equation
81  template<class Type>
82  void addSupType
83  (
84  const VolField<Type>& field,
85  fvMatrix<Type>& eqn
86  ) const;
87 
88  //- Add a source term to a compressible equation
89  template<class Type>
90  void addSupType
91  (
92  const volScalarField& rho,
93  const VolField<Type>& field,
94  fvMatrix<Type>& eqn
95  ) const;
96 
97  //- Add a source term to a phase equation
98  template<class Type>
99  void addSupType
100  (
101  const volScalarField& alpha,
102  const volScalarField& rho,
103  const VolField<Type>& field,
104  fvMatrix<Type>& eqn
105  ) const;
106 
107 
108 public:
109 
110  //- Runtime type information
111  TypeName("fvTotalSource");
112 
113 
114  // Constructors
115 
116  //- Construct from explicit source name and mesh
118  (
119  const word& name,
120  const word& modelType,
121  const fvMesh& mesh,
122  const dictionary& dict
123  );
124 
125 
126  //- Destructor
127  virtual ~fvTotalSource();
128 
129 
130  // Member Functions
131 
132  // Checks
133 
134  //- Return true if the fvModel adds a source term to the given
135  // field's transport equation
136  virtual bool addsSupToField(const word& fieldName) const;
137 
138 
139  // Access
140 
141  //- Return the phase name
142  inline const word& phaseName() const;
143 
144  //- Return the volume of cells that the source applies to
145  virtual scalar V() const = 0;
146 
147 
148  // Sources
149 
150  //- Return the source value
151  virtual dimensionedScalar S() const = 0;
152 
153  //- Return the source value
154  virtual tmp<scalarField> source(const word& fieldName) const;
155 
156 
157  // IO
158 
159  //- Read source dictionary
160  virtual bool read(const dictionary& dict);
161 };
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #include "fvTotalSourceI.H"
171 
172 #ifdef NoRepository
173  #include "fvTotalSourceTemplates.C"
174 #endif
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 #endif
179 
180 // ************************************************************************* //
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:47
static const word & fieldName()
Return the name of the field associated with a source term. Special.
Definition: fvModelI.H:39
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:55
void addSupType(const VolField< Type > &field, fvMatrix< Type > &eqn) const
Add a source term to an equation.
virtual scalar V() const =0
Return the volume of cells that the source applies to.
const word & phaseName() const
Return the phase name.
TypeName("fvTotalSource")
Runtime type information.
fvTotalSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: fvTotalSource.C:73
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual tmp< scalarField > source(const word &fieldName) const
Return the source value.
virtual ~fvTotalSource()
Destructor.
Definition: fvTotalSource.C:89
virtual dimensionedScalar S() const =0
Return the source value.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: fvTotalSource.C:95
void addSource(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: fvTotalSource.C:48
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
dictionary dict