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-2021 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 volume normalisation
68  switch (volumeMode_)
69  {
71  VDash_ = set_.V();
72  break;
73  case volumeMode::specific:
74  VDash_ = 1;
75  break;
76  }
77 
78  // Set field source terms
79  fieldSp_.clear();
80  fieldSu_.clear();
81  forAllConstIter(dictionary, coeffs().subDict("sources"), iter)
82  {
83  fieldSu_.set
84  (
85  iter().keyword(),
86  objectFunction1::New<VolField>
87  (
88  "explicit",
89  iter().dict(),
90  iter().keyword(),
91  mesh(),
92  false
93  ).ptr()
94  );
95  fieldSp_.set
96  (
97  iter().keyword(),
99  (
100  "implicit",
101  iter().dict()
102  ).ptr()
103  );
104  }
105 }
106 
107 
108 template<class Type>
109 void Foam::fv::semiImplicitSource::addSupType
110 (
111  fvMatrix<Type>& eqn,
112  const word& fieldName
113 ) const
114 {
115  if (debug)
116  {
117  Info<< "semiImplicitSource<" << pTraits<Type>::typeName
118  << ">::addSup for source " << name() << endl;
119  }
120 
121  const scalar t = mesh().time().value();
122 
124 
126  (
127  IOobject
128  (
129  name() + fieldName + "Su",
130  mesh().time().timeName(),
131  mesh(),
134  ),
135  mesh(),
137  (
138  "zero",
139  eqn.dimensions()/dimVolume,
140  Zero
141  ),
142  false
143  );
144 
145  // Explicit source function for the field
146  UIndirectList<Type>(Su, set_.cells()) =
147  fieldSu_[fieldName]->value<Type>(t)/VDash_;
148 
150  (
151  IOobject
152  (
153  name() + fieldName + "Sp",
154  mesh().time().timeName(),
155  mesh(),
158  ),
159  mesh(),
161  (
162  "zero",
163  Su.dimensions()/psi.dimensions(),
164  0
165  ),
166  false
167  );
168 
169  // Implicit source function for the field
170  UIndirectList<scalar>(Sp, set_.cells()) =
171  fieldSp_[fieldName]->value(t)/VDash_;
172 
173  eqn += Su + fvm::SuSp(Sp, psi);
174 }
175 
176 
177 template<class Type>
178 void Foam::fv::semiImplicitSource::addSupType
179 (
180  const volScalarField& rho,
181  fvMatrix<Type>& eqn,
182  const word& fieldName
183 ) const
184 {
185  return this->addSup(eqn, fieldName);
186 }
187 
188 
189 template<class Type>
190 void Foam::fv::semiImplicitSource::addSupType
191 (
192  const volScalarField& alpha,
193  const volScalarField& rho,
194  fvMatrix<Type>& eqn,
195  const word& fieldName
196 ) const
197 {
198  return this->addSup(eqn, fieldName);
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
203 
205 (
206  const word& name,
207  const word& modelType,
208  const dictionary& dict,
209  const fvMesh& mesh
210 )
211 :
212  fvModel(name, modelType, dict, mesh),
213  set_(coeffs(), mesh),
214  volumeMode_(volumeMode::absolute),
215  VDash_(1)
216 {
217  readCoeffs();
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
222 
224 {}
225 
226 
227 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 
230 {
231  return fieldSu_.toc();
232 }
233 
234 
236 
237 
239 
240 
242 (
245 );
246 
247 
249 {
250  set_.updateMesh(mpm);
251 }
252 
253 
255 {
256  if (fvModel::read(dict))
257  {
258  set_.read(coeffs());
259  readCoeffs();
260  return true;
261  }
262  else
263  {
264  return false;
265  }
266 }
267 
268 
269 // ************************************************************************* //
#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)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:158
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:55
Generic GeometricField class.
Generic dimensioned Type class.
Semi-implicit source, described using an input dictionary. The injection rate coefficients are specif...
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.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
FOR_ALL_FIELD_TYPES(IMPLEMENT_FV_MODEL_ADD_SUP, fv::semiImplicitSource)
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
virtual bool read(const dictionary &dict)
Read source dictionary.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
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
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
#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:188
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:78
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.
virtual void updateMesh(const mapPolyMesh &)
Update for mesh changes.
const dimensionSet dimVolume
messageStream Info
const volScalarField & psi
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
virtual ~semiImplicitSource()
Destructor.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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.
const dimensionSet & dimensions() const
Definition: fvMatrix.H:287