semiImplicitSource.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) 2011-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::fv::semiImplicitSource
26 
27 Description
28  Semi-implicit source, described using an input dictionary. The injection
29  rate coefficients are specified as pairs of Su-Sp coefficients, i.e.
30 
31  \f[
32  S(x) = S_u + S_p x
33  \f]
34 
35  where
36  \vartable
37  S(x) | net source for field 'x'
38  S_u | explicit source contribution
39  S_p | linearised implicit contribution
40  \endvartable
41 
42  Example tabulated heat source specification for internal energy:
43  \verbatim
44  volumeMode absolute; // specific
45  sources
46  {
47  e
48  {
49  explicit table ((0 0) (1.5 $power));
50  implicit 0;
51  }
52  }
53  \endverbatim
54 
55  Example coded heat source specification for enthalpy:
56  \verbatim
57  volumeMode absolute; // specific
58  sources
59  {
60  h
61  {
62  explicit
63  {
64  type coded;
65  name heatInjection;
66  code
67  #{
68  // Power amplitude
69  const scalar powerAmplitude = 1000;
70 
71  // x is the current time
72  return mag(powerAmplitude*sin(x));
73  #};
74  }
75  implicit 0;
76  }
77  }
78  \endverbatim
79 
80  Valid fvModels for the \c volumeMode entry include:
81  - absolute: values are given as <quantity>
82  - specific: values are given as <quantity>/m3
83 
84 See also
85  Foam::fvModel
86 
87 SourceFiles
88  semiImplicitSource.C
89 
90 \*---------------------------------------------------------------------------*/
91 
92 #ifndef semiImplicitSource_H
93 #define semiImplicitSource_H
94 
95 #include "fvModel.H"
96 #include "fvCellSet.H"
97 #include "HashPtrTable.H"
98 #include "unknownTypeFunction1.H"
99 
100 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102 namespace Foam
103 {
104 namespace fv
105 {
106 
107 /*---------------------------------------------------------------------------*\
108  Class semiImplicitSource Declaration
109 \*---------------------------------------------------------------------------*/
110 
111 class semiImplicitSource
112 :
113  public fvModel
114 {
115 public:
116 
117  // Public data
118 
119  //- Enumeration for volume types
120  enum class volumeMode
121  {
123  specific
124  };
125 
126  //- Property type names
128 
129 
130 private:
131 
132  // Private member data
133 
134  //- The set of cells the fvConstraint applies to
135  fvCellSet set_;
136 
137  //- Volume mode
138  volumeMode volumeMode_;
139 
140  //- Explicit parts of the sources
142 
143  //- Implicit parts of the sources
145 
146 
147  // Private Member Functions
148 
149  //- Non-virtual read
150  void readCoeffs();
151 
152 
153  // Sources
154 
155  //- Add a source term to an equation
156  template <class Type>
157  void addSupType(fvMatrix<Type>& eqn, const word& fieldName) const;
158 
159  //- Add a source term to a compressible equation
160  template <class Type>
161  void addSupType
162  (
163  const volScalarField& rho,
164  fvMatrix<Type>& eqn,
165  const word& fieldName
166  ) const;
167 
168  //- Add a source term to a phase equation
169  template <class Type>
170  void addSupType
171  (
172  const volScalarField& alpha,
173  const volScalarField& rho,
174  fvMatrix<Type>& eqn,
175  const word& fieldName
176  ) const;
177 
178 
179 public:
180 
181  //- Runtime type information
182  TypeName("semiImplicitSource");
183 
184 
185  // Constructors
186 
187  //- Construct from components
189  (
190  const word& name,
191  const word& modelType,
192  const fvMesh& mesh,
194  );
195 
196  //- Destructor
197  virtual ~semiImplicitSource();
198 
199 
200  // Member Functions
201 
202  // Checks
203 
204  //- Return the list of fields for which the fvModel adds source term
205  // to the transport equation
206  virtual wordList addSupFields() const;
207 
208 
209  // Sources
210 
211  //- Add a source term to an equation
213 
214  //- Add a source term to a compressible equation
216 
217  //- Add a source term to a phase equation
219 
220 
221  // Mesh changes
222 
223  //- Update for mesh motion
224  virtual bool movePoints();
225 
226  //- Update topology using the given map
227  virtual void topoChange(const polyTopoChangeMap&);
228 
229  //- Update from another mesh using the given map
230  virtual void mapMesh(const polyMeshMap&);
231 
232  //- Redistribute or update using the given distribution map
233  virtual void distribute(const polyDistributionMap&);
234 
235 
236  // IO
237 
238  //- Read source dictionary
239  virtual bool read(const dictionary& dict);
240 };
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 } // End namespace fv
246 } // End namespace Foam
247 
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
Generic GeometricField class.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
TypeName("semiImplicitSource")
Runtime type information.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual ~semiImplicitSource()
Destructor.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
volumeMode
Enumeration for volume types.
FOR_ALL_FIELD_TYPES(DEFINE_FV_MODEL_ADD_SUP)
Add a source term to an equation.
semiImplicitSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
static const NamedEnum< volumeMode, 2 > volumeModeNames_
Property type names.
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
#define DEFINE_FV_MODEL_ADD_SUP(Type, nullArg)
Definition: fvModelM.H:26
#define DEFINE_FV_MODEL_ADD_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:43
#define DEFINE_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, nullArg)
Definition: fvModelM.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
labelList fv(nPoints)
dictionary dict