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-2020 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  cellSetOption,
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 template<class Type>
63 void Foam::fv::semiImplicitSource::addSupType
64 (
65  fvMatrix<Type>& eqn,
66  const label fieldi
67 )
68 {
69  if (debug)
70  {
71  Info<< "semiImplicitSource<" << pTraits<Type>::typeName
72  << ">::addSup for source " << name_ << endl;
73  }
74 
75  const scalar t = mesh_.time().value();
76 
78 
80  (
81  IOobject
82  (
83  name_ + fieldNames_[fieldi] + "Su",
84  mesh_.time().timeName(),
85  mesh_,
88  ),
89  mesh_,
91  (
92  "zero",
93  eqn.dimensions()/dimVolume,
94  Zero
95  ),
96  false
97  );
98 
99  // Explicit source function for the field
100  UIndirectList<Type>(Su, cells_) = fieldSu_[fieldi].value<Type>(t)/VDash_;
101 
103  (
104  IOobject
105  (
106  name_ + fieldNames_[fieldi] + "Sp",
107  mesh_.time().timeName(),
108  mesh_,
111  ),
112  mesh_,
114  (
115  "zero",
116  Su.dimensions()/psi.dimensions(),
117  0
118  ),
119  false
120  );
121 
122  // Implicit source function for the field
123  UIndirectList<scalar>(Sp, cells_) = fieldSp_[fieldi].value(t)/VDash_;
124 
125  eqn += Su + fvm::SuSp(Sp, psi);
126 }
127 
128 
129 template<class Type>
130 void Foam::fv::semiImplicitSource::addSupType
131 (
132  const volScalarField& rho,
133  fvMatrix<Type>& eqn,
134  const label fieldi
135 )
136 {
137  if (debug)
138  {
139  Info<< "semiImplicitSource<" << pTraits<Type>::typeName
140  << ">::addSup for source " << name_ << endl;
141  }
142 
143  return this->addSup(eqn, fieldi);
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 
150 (
151  const word& name,
152  const word& modelType,
153  const dictionary& dict,
154  const fvMesh& mesh
155 )
156 :
157  cellSetOption(name, modelType, dict, mesh),
158  volumeMode_(volumeMode::absolute),
159  VDash_(1)
160 {
161  read(dict);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
166 
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 (
175  fvMatrix<scalar>& eqn,
176  const label fieldi
177 )
178 {
179  addSupType(eqn, fieldi);
180 }
181 
182 
184 (
185  fvMatrix<vector>& eqn,
186  const label fieldi
187 )
188 {
189  addSupType(eqn, fieldi);
190 }
191 
192 
194 (
196  const label fieldi
197 )
198 {
199  addSupType(eqn, fieldi);
200 }
201 
202 
204 (
206  const label fieldi
207 )
208 {
209  addSupType(eqn, fieldi);
210 }
211 
212 
214 (
215  fvMatrix<tensor>& eqn,
216  const label fieldi
217 )
218 {
219  addSupType(eqn, fieldi);
220 }
221 
222 
224 (
225  const volScalarField& rho,
226  fvMatrix<scalar>& eqn,
227  const label fieldi
228 )
229 {
230  addSupType(eqn, fieldi);
231 }
232 
233 
235 (
236  const volScalarField& rho,
237  fvMatrix<vector>& eqn,
238  const label fieldi
239 )
240 {
241  addSupType(eqn, fieldi);
242 }
243 
244 
246 (
247  const volScalarField& rho,
249  const label fieldi
250 )
251 {
252  addSupType(eqn, fieldi);
253 }
254 
255 
257 (
258  const volScalarField& rho,
260  const label fieldi
261 )
262 {
263  addSupType(eqn, fieldi);
264 }
265 
266 
268 (
269  const volScalarField& rho,
270  fvMatrix<tensor>& eqn,
271  const label fieldi
272 )
273 {
274  addSupType(eqn, fieldi);
275 }
276 
277 
279 (
280  const volScalarField& alpha,
281  const volScalarField& rho,
282  fvMatrix<scalar>& eqn,
283  const label fieldi
284 )
285 {
286  addSupType(eqn, fieldi);
287 }
288 
289 
291 (
292  const volScalarField& alpha,
293  const volScalarField& rho,
294  fvMatrix<vector>& eqn,
295  const label fieldi
296 )
297 {
298  addSupType(eqn, fieldi);
299 }
300 
301 
303 (
304  const volScalarField& alpha,
305  const volScalarField& rho,
307  const label fieldi
308 )
309 {
310  addSupType(eqn, fieldi);
311 }
312 
313 
315 (
316  const volScalarField& alpha,
317  const volScalarField& rho,
319  const label fieldi
320 )
321 {
322  addSupType(eqn, fieldi);
323 }
324 
325 
327 (
328  const volScalarField& alpha,
329  const volScalarField& rho,
330  fvMatrix<tensor>& eqn,
331  const label fieldi
332 )
333 {
334  addSupType(eqn, fieldi);
335 }
336 
337 
339 {
340  if (cellSetOption::read(dict))
341  {
342  volumeMode_ = volumeModeNames_.read(coeffs_.lookup("volumeMode"));
343 
344  const dictionary& sources = coeffs_.subDict("sources");
345 
346  // Number of fields with a source term
347  const label nFields = sources.size();
348 
349  // Set field names and source terms
350  fieldNames_.setSize(nFields);
351  fieldSp_.setSize(nFields);
352  fieldSu_.setSize(nFields);
353  label i = 0;
354  forAllConstIter(dictionary, sources, iter)
355  {
356  fieldNames_[i] = iter().keyword();
357  fieldSu_.set
358  (
359  i,
360  objectFunction1::New<VolField>
361  (
362  "explicit",
363  iter().dict(),
364  fieldNames_[i],
365  mesh_
366  ).ptr()
367  );
368  fieldSp_.set
369  (
370  i,
372  (
373  "implicit",
374  iter().dict()
375  ).ptr()
376  );
377  i++;
378  }
379 
380  // Set volume normalisation
381  if (volumeMode_ == volumeMode::absolute)
382  {
383  VDash_ = V_;
384  }
385 
386  applied_.setSize(nFields, false);
387 
388  return true;
389  }
390  else
391  {
392  return false;
393  }
394 }
395 
396 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:52
dictionary dict
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
zeroField Su
Definition: alphaSuSp.H:1
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Traits class for primitives.
Definition: pTraits.H:50
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
virtual bool read(const dictionary &dict)
Read source dictionary.
Generic GeometricField class.
Generic dimensioned Type class.
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 dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
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
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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.
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.
messageStream Info
Cell-set options abstract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:69
const volScalarField & psi
virtual ~semiImplicitSource()
Destructor.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
zeroField Sp
Definition: alphaSuSp.H:2
const dimensionSet & dimensions() const
Definition: fvMatrix.H:287