SemiImplicitSource.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "SemiImplicitSource.H"
29 #include "fvMesh.H"
30 #include "fvMatrices.H"
31 #include "fvmSup.H"
32 
33 // * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
37 {
38  "absolute", "specific"
39 };
40 
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
44 template<class Type>
47 (
48  const word& vmtName
49 ) const
50 {
51  forAll(volumeModeTypeNames_, i)
52  {
53  if (vmtName == volumeModeTypeNames_[i])
54  {
55  return volumeModeType(i);
56  }
57  }
58 
60  << "Unknown volumeMode type " << vmtName
61  << ". Valid volumeMode types are:" << nl << volumeModeTypeNames_
62  << exit(FatalError);
63 
64  return volumeModeType(0);
65 }
66 
67 
68 template<class Type>
70 (
71  const volumeModeType& vmtType
72 ) const
73 {
74  if (vmtType > volumeModeTypeNames_.size())
75  {
76  return "UNKNOWN";
77  }
78  else
79  {
80  return volumeModeTypeNames_[vmtType];
81  }
82 }
83 
84 
85 template<class Type>
87 {
88  label count = dict.size();
89 
90  fieldNames_.resize(count);
91  injectionRate_.resize(count);
92  applied_.resize(count, false);
93 
94  count = 0;
95  for (const entry& dEntry : dict)
96  {
97  fieldNames_[count] = dEntry.keyword();
98  dEntry.readEntry(injectionRate_[count]);
99 
100  ++count;
101  }
102 
103  // Set volume normalisation
104  if (volumeMode_ == vmAbsolute)
105  {
106  VDash_ = V_;
107  }
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
113 template<class Type>
115 (
116  const word& name,
117  const word& modelType,
118  const dictionary& dict,
119  const fvMesh& mesh
120 )
121 :
122  cellSetOption(name, modelType, dict, mesh),
123  volumeMode_(vmAbsolute),
124  VDash_(1.0),
125  injectionRate_()
126 {
127  read(dict);
128 }
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
133 template<class Type>
135 (
136  fvMatrix<Type>& eqn,
137  const label fieldi
138 )
139 {
140  if (debug)
141  {
142  Info<< "SemiImplicitSource<" << pTraits<Type>::typeName
143  << ">::addSup for source " << name_ << endl;
144  }
145 
147 
149  (
150  IOobject
151  (
152  name_ + fieldNames_[fieldi] + "Su",
153  mesh_.time().timeName(),
154  mesh_,
155  IOobject::NO_READ,
156  IOobject::NO_WRITE
157  ),
158  mesh_,
160  false
161  );
162 
163  UIndirectList<Type>(Su, cells_) = injectionRate_[fieldi].first()/VDash_;
164 
166  (
167  IOobject
168  (
169  name_ + fieldNames_[fieldi] + "Sp",
170  mesh_.time().timeName(),
171  mesh_,
172  IOobject::NO_READ,
173  IOobject::NO_WRITE
174  ),
175  mesh_,
176  dimensioned<scalar>(Su.dimensions()/psi.dimensions(), Zero),
177  false
178  );
179 
180  UIndirectList<scalar>(Sp, cells_) = injectionRate_[fieldi].second()/VDash_;
181 
182  eqn += Su + fvm::SuSp(Sp, psi);
183 }
184 
185 
186 template<class Type>
188 (
189  const volScalarField& rho,
190  fvMatrix<Type>& eqn,
191  const label fieldi
192 )
193 {
194  if (debug)
195  {
196  Info<< "SemiImplicitSource<" << pTraits<Type>::typeName
197  << ">::addSup for source " << name_ << endl;
198  }
199 
200  return this->addSup(eqn, fieldi);
201 }
202 
203 
204 
205 template<class Type>
207 {
209  {
210  volumeMode_ = wordToVolumeModeType(coeffs_.get<word>("volumeMode"));
211  setFieldData(coeffs_.subDict("injectionRateSuSp"));
212 
213  return true;
214  }
215 
216  return false;
217 }
218 
219 
220 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::cellSetOption
Cell-set options abtract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:72
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::fvMatrix::dimensions
const dimensionSet & dimensions() const
Definition: fvMatrix.H:290
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::fv::SemiImplicitSource::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: SemiImplicitSource.C:206
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
rho
rho
Definition: readInitialConditions.H:96
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:285
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::fam::SuSp
tmp< faMatrix< Type > > SuSp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famSup.C:151
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::fv::SemiImplicitSource
Semi-implicit source, described using an input dictionary. The injection rate coefficients are specif...
Definition: SemiImplicitSource.H:95
Foam::label
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:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fv::SemiImplicitSource::volumeModeTypeToWord
word volumeModeTypeToWord(const volumeModeType &vtType) const
Helper function to convert from a volumeModeType to a word.
Definition: SemiImplicitSource.C:70
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::fv::SemiImplicitSource::setFieldData
void setFieldData(const dictionary &dict)
Set the local field data.
Definition: SemiImplicitSource.C:86
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:74
Foam::List< word >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
SemiImplicitSource.H
Foam::fv::SemiImplicitSource::addSup
virtual void addSup(fvMatrix< Type > &eqn, const label fieldi)
Add explicit contribution to equation.
Definition: SemiImplicitSource.C:135
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
Foam::fv::SemiImplicitSource::SemiImplicitSource
SemiImplicitSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: SemiImplicitSource.C:115
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< Type, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fv::SemiImplicitSource::volumeModeType
volumeModeType
Enumeration for volume types.
Definition: SemiImplicitSource.H:119
Foam::fv::SemiImplicitSource::wordToVolumeModeType
volumeModeType wordToVolumeModeType(const word &vtName) const
Helper function to convert from a word to a volumeModeType.
Definition: SemiImplicitSource.C:47