semiImplicitSource.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) 2011-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 "semiImplicitSource.H"
27 #include "fvMesh.H"
28 #include "fvMatrices.H"
29 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
39 
41  (
42  fvModel,
45  );
46  }
47 }
48 
51 {
52  "absolute",
53  "specific"
54 };
55 
56 
57 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58 
59 void Foam::fv::semiImplicitSource::readCoeffs(const dictionary& dict)
60 {
61  // Get the volume mode
62  volumeMode_ = volumeModeNames_.read(dict.lookup("volumeMode"));
63 
64  // Set field source terms
65  fieldSu_.clear();
66  fieldSp_.clear();
67  forAllConstIter(dictionary, dict.subDict("sources"), iter)
68  {
69  fieldSu_.set
70  (
71  iter().keyword(),
72  new unknownTypeFunction1
73  (
74  "explicit",
75  mesh().time().userUnits(),
76  iter().dict()
77  )
78  );
79  fieldSp_.set
80  (
81  iter().keyword(),
82  new unknownTypeFunction1
83  (
84  "implicit",
85  mesh().time().userUnits(),
86  iter().dict()
87  )
88  );
89  }
90 }
91 
92 
93 template<class Type>
94 void Foam::fv::semiImplicitSource::addSupType
95 (
96  const VolField<Type>& field,
97  fvMatrix<Type>& eqn
98 ) const
99 {
100  // Set the value units for the functions
101  fieldSu_[field.name()]->template setValueUnits<Type>
102  (
103  eqn.dimensions()
104  );
105  fieldSp_[field.name()]->template setValueUnits<scalar>
106  (
107  eqn.dimensions()/eqn.psi().dimensions()
108  );
109 
110  const scalar t = mesh().time().value();
111 
112  const VolField<Type>& psi = eqn.psi();
113 
114  VolInternalField<Type> Su
115  (
116  IOobject
117  (
118  name() + field.name() + "Su",
119  mesh().time().name(),
120  mesh(),
123  ),
124  mesh(),
125  dimensioned<Type>
126  (
127  "zero",
128  eqn.dimensions()/dimVolume,
129  Zero
130  ),
131  false
132  );
133 
134  // Set volume normalisation
135  scalar VDash = NaN;
136  switch (volumeMode_)
137  {
139  VDash = zone_.V();
140  break;
141  case volumeMode::specific:
142  VDash = 1;
143  break;
144  }
145 
146  // Explicit source function for the field
147  UIndirectList<Type>(Su, zone_.zone()) =
148  fieldSu_[field.name()]->template value<Type>(t)/VDash;
149 
151  (
152  IOobject
153  (
154  name() + field.name() + "Sp",
155  mesh().time().name(),
156  mesh(),
159  ),
160  mesh(),
161  dimensioned<scalar>
162  (
163  "zero",
164  Su.dimensions()/psi.dimensions(),
165  0
166  ),
167  false
168  );
169 
170  // Implicit source function for the field
171  UIndirectList<scalar>(Sp, zone_.zone()) =
172  fieldSp_[field.name()]->template value<scalar>(t)/VDash;
173 
174  eqn += Su - fvm::SuSp(-Sp, psi);
175 }
176 
177 
178 template<class Type>
179 void Foam::fv::semiImplicitSource::addSupType
180 (
181  const volScalarField& rho,
182  const VolField<Type>& field,
183  fvMatrix<Type>& eqn
184 ) const
185 {
186  return addSup(field, eqn);
187 }
188 
189 
190 template<class Type>
191 void Foam::fv::semiImplicitSource::addSupType
192 (
193  const volScalarField& alpha,
194  const volScalarField& rho,
195  const VolField<Type>& field,
196  fvMatrix<Type>& eqn
197 ) const
198 {
199  return addSup(field, eqn);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
204 
206 (
207  const word& name,
208  const word& modelType,
209  const fvMesh& mesh,
210  const dictionary& dict
211 )
212 :
213  fvModel(name, modelType, mesh, dict),
214  zone_(mesh, coeffs(dict)),
215  volumeMode_(volumeMode::absolute)
216 {
217  readCoeffs(coeffs(dict));
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
222 
224 {}
225 
226 
227 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 
230 {
231  return fieldSu_.toc();
232 }
233 
234 
236 (
238  fv::semiImplicitSource
239 )
240 
241 
243 (
245  fv::semiImplicitSource
246 )
247 
248 
250 (
252  fv::semiImplicitSource
253 )
254 
255 
257 {
258  zone_.movePoints();
259  return true;
260 }
261 
262 
264 {
265  zone_.topoChange(map);
266 }
267 
268 
270 {
271  zone_.mapMesh(map);
272 }
273 
274 
276 (
277  const polyDistributionMap& map
278 )
279 {
280  zone_.distribute(map);
281 }
282 
283 
285 {
286  if (fvModel::read(dict))
287  {
288  zone_.read(coeffs(dict));
289  readCoeffs(coeffs(dict));
290  return true;
291  }
292  else
293  {
294  return false;
295  }
296 }
297 
298 
299 // ************************************************************************* //
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:54
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Finite volume model abstract base class.
Definition: fvModel.H:60
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
const word & keyword() const
Return name as the keyword.
Definition: fvModelI.H:63
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
Semi-implicit source, described using an input dictionary. The injection rate coefficients are specif...
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
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 movePoints()
Add a source term to an equation.
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.
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
Calculate the matrix for implicit and explicit sources.
const volScalarField & psi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
labelList fv(nPoints)
dictionary dict