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-2022 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  {
38  defineTypeNameAndDebug(semiImplicitSource, 0);
39 
41  (
42  fvModel,
43  semiImplicitSource,
44  dictionary
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("explicit", iter().dict())
76  );
77  fieldSp_.set
78  (
79  iter().keyword(),
80  Function1<scalar>::New("implicit", iter().dict()).ptr()
81  );
82  }
83 }
84 
85 
86 template<class Type>
87 void Foam::fv::semiImplicitSource::addSupType
88 (
89  fvMatrix<Type>& eqn,
90  const word& fieldName
91 ) const
92 {
93  if (debug)
94  {
95  Info<< "semiImplicitSource<" << pTraits<Type>::typeName
96  << ">::addSup for source " << name() << endl;
97  }
98 
99  const scalar t = mesh().time().userTimeValue();
100 
102 
104  (
105  IOobject
106  (
107  name() + fieldName + "Su",
108  mesh().time().timeName(),
109  mesh(),
112  ),
113  mesh(),
115  (
116  "zero",
117  eqn.dimensions()/dimVolume,
118  Zero
119  ),
120  false
121  );
122 
123  // Set volume normalisation
124  scalar VDash = NaN;
125  switch (volumeMode_)
126  {
128  VDash = set_.V();
129  break;
130  case volumeMode::specific:
131  VDash = 1;
132  break;
133  }
134 
135  // Explicit source function for the field
136  UIndirectList<Type>(Su, set_.cells()) =
137  fieldSu_[fieldName]->value<Type>(t)/VDash;
138 
140  (
141  IOobject
142  (
143  name() + fieldName + "Sp",
144  mesh().time().timeName(),
145  mesh(),
148  ),
149  mesh(),
151  (
152  "zero",
153  Su.dimensions()/psi.dimensions(),
154  0
155  ),
156  false
157  );
158 
159  // Implicit source function for the field
160  UIndirectList<scalar>(Sp, set_.cells()) =
161  fieldSp_[fieldName]->value(t)/VDash;
162 
163  eqn += Su + fvm::SuSp(Sp, psi);
164 }
165 
166 
167 template<class Type>
168 void Foam::fv::semiImplicitSource::addSupType
169 (
170  const volScalarField& rho,
171  fvMatrix<Type>& eqn,
172  const word& fieldName
173 ) const
174 {
175  return this->addSup(eqn, fieldName);
176 }
177 
178 
179 template<class Type>
180 void Foam::fv::semiImplicitSource::addSupType
181 (
182  const volScalarField& alpha,
183  const volScalarField& rho,
184  fvMatrix<Type>& eqn,
185  const word& fieldName
186 ) const
187 {
188  return this->addSup(eqn, fieldName);
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
193 
195 (
196  const word& name,
197  const word& modelType,
198  const dictionary& dict,
199  const fvMesh& mesh
200 )
201 :
202  fvModel(name, modelType, dict, mesh),
203  set_(coeffs(), mesh),
204  volumeMode_(volumeMode::absolute)
205 {
206  readCoeffs();
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
211 
213 {}
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
219 {
220  return fieldSu_.toc();
221 }
222 
223 
225 
226 
228 
229 
231 (
234 );
235 
236 
238 {
239  set_.movePoints();
240  return true;
241 }
242 
243 
245 {
246  set_.topoChange(map);
247 }
248 
249 
251 {
252  set_.mapMesh(map);
253 }
254 
255 
257 (
258  const polyDistributionMap& map
259 )
260 {
261  set_.distribute(map);
262 }
263 
264 
266 {
267  if (fvModel::read(dict))
268  {
269  set_.read(coeffs());
270  readCoeffs();
271  return true;
272  }
273  else
274  {
275  return false;
276  }
277 }
278 
279 
280 // ************************************************************************* //
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
Run-time selectable general function of one variable.
Definition: Function1.H:52
dictionary dict
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Wrapper around Function1 that constructs a function for an as yet unknown primitive type...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:290
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Traits class for primitives.
Definition: pTraits.H:50
Finite volume model abstract base class.
Definition: fvModel.H:57
Generic GeometricField class.
Generic dimensioned Type class.
Semi-implicit source, described using an input dictionary. The injection rate coefficients are specif...
virtual bool movePoints()
Update for mesh motion.
fvMesh & mesh
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
semiImplicitSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual bool read(const dictionary &dict)
Read source dictionary.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const volScalarField & psi
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
Definition: fvModelM.H:33
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
static const NamedEnum< volumeMode, 2 > volumeModeNames_
Property type names.
#define IMPLEMENT_FV_MODEL_ADD_RHO_SUP(Type, modelType)
Definition: fvModelM.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
A List with indirect addressing.
Definition: fvMatrix.H:106
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const dimensionSet dimVolume
messageStream Info
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
virtual ~semiImplicitSource()
Destructor.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
const dimensionSet & dimensions() const
Definition: fvMatrix.H:295