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-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 "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  template<>
50  {
51  "absolute",
52  "specific"
53  };
54 }
55 
58 
59 
60 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61 
62 void Foam::fv::semiImplicitSource::readCoeffs()
63 {
64  // Get the volume mode
65  volumeMode_ = volumeModeNames_.read(coeffs().lookup("volumeMode"));
66 
67  // Set field source terms
68  fieldSu_.clear();
69  fieldSp_.clear();
70  forAllConstIter(dictionary, coeffs().subDict("sources"), iter)
71  {
72  fieldSu_.set
73  (
74  iter().keyword(),
75  new unknownTypeFunction1
76  (
77  "explicit",
78  mesh().time().userUnits(),
79  iter().dict()
80  )
81  );
82  fieldSp_.set
83  (
84  iter().keyword(),
85  new unknownTypeFunction1
86  (
87  "implicit",
88  mesh().time().userUnits(),
89  iter().dict()
90  )
91  );
92  }
93 }
94 
95 
96 template<class Type>
97 void Foam::fv::semiImplicitSource::addSupType
98 (
99  const VolField<Type>& field,
100  fvMatrix<Type>& eqn
101 ) const
102 {
103  // Set the value units for the functions
104  fieldSu_[field.name()]->template setValueUnits<Type>
105  (
106  eqn.dimensions()
107  );
108  fieldSp_[field.name()]->template setValueUnits<scalar>
109  (
110  eqn.dimensions()/eqn.psi().dimensions()
111  );
112 
113  const scalar t = mesh().time().value();
114 
115  const VolField<Type>& psi = eqn.psi();
116 
117  VolInternalField<Type> Su
118  (
119  IOobject
120  (
121  name() + field.name() + "Su",
122  mesh().time().name(),
123  mesh(),
126  ),
127  mesh(),
128  dimensioned<Type>
129  (
130  "zero",
131  eqn.dimensions()/dimVolume,
132  Zero
133  ),
134  false
135  );
136 
137  // Set volume normalisation
138  scalar VDash = NaN;
139  switch (volumeMode_)
140  {
142  VDash = set_.V();
143  break;
144  case volumeMode::specific:
145  VDash = 1;
146  break;
147  }
148 
149  // Explicit source function for the field
150  UIndirectList<Type>(Su, set_.cells()) =
151  fieldSu_[field.name()]->template value<Type>(t)/VDash;
152 
154  (
155  IOobject
156  (
157  name() + field.name() + "Sp",
158  mesh().time().name(),
159  mesh(),
162  ),
163  mesh(),
164  dimensioned<scalar>
165  (
166  "zero",
167  Su.dimensions()/psi.dimensions(),
168  0
169  ),
170  false
171  );
172 
173  // Implicit source function for the field
174  UIndirectList<scalar>(Sp, set_.cells()) =
175  fieldSp_[field.name()]->template value<scalar>(t)/VDash;
176 
177  eqn += Su - fvm::SuSp(-Sp, psi);
178 }
179 
180 
181 template<class Type>
182 void Foam::fv::semiImplicitSource::addSupType
183 (
184  const volScalarField& rho,
185  const VolField<Type>& field,
186  fvMatrix<Type>& eqn
187 ) const
188 {
189  return addSup(field, eqn);
190 }
191 
192 
193 template<class Type>
194 void Foam::fv::semiImplicitSource::addSupType
195 (
196  const volScalarField& alpha,
197  const volScalarField& rho,
198  const VolField<Type>& field,
199  fvMatrix<Type>& eqn
200 ) const
201 {
202  return addSup(field, eqn);
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
207 
209 (
210  const word& name,
211  const word& modelType,
212  const fvMesh& mesh,
213  const dictionary& dict
214 )
215 :
216  fvModel(name, modelType, mesh, dict),
217  set_(mesh, coeffs()),
218  volumeMode_(volumeMode::absolute)
219 {
220  readCoeffs();
221 }
222 
223 
224 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
225 
227 {}
228 
229 
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
231 
233 {
234  return fieldSu_.toc();
235 }
236 
237 
239 (
241  fv::semiImplicitSource
242 )
243 
244 
246 (
248  fv::semiImplicitSource
249 )
250 
251 
253 (
255  fv::semiImplicitSource
256 )
257 
258 
260 {
261  set_.movePoints();
262  return true;
263 }
264 
265 
267 {
268  set_.topoChange(map);
269 }
270 
271 
273 {
274  set_.mapMesh(map);
275 }
276 
277 
279 (
280  const polyDistributionMap& map
281 )
282 {
283  set_.distribute(map);
284 }
285 
286 
288 {
289  if (fvModel::read(dict))
290  {
291  set_.read(coeffs());
292  readCoeffs();
293  return true;
294  }
295  else
296  {
297  return false;
298  }
299 }
300 
301 
302 // ************************************************************************* //
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
DimensionedField< Type, GeoMesh > 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:54
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
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
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
labelList fv(nPoints)
dictionary dict