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 | Copyright (C) 2011-2016 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"
30 
31 // * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  IStringStream("(absolute specific)")()
37 );
38 
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
42 template<class Type>
45 (
46  const word& vmtName
47 ) const
48 {
49  forAll(volumeModeTypeNames_, i)
50  {
51  if (vmtName == volumeModeTypeNames_[i])
52  {
53  return volumeModeType(i);
54  }
55  }
56 
58  << "Unknown volumeMode type " << vmtName
59  << ". Valid volumeMode types are:" << nl << volumeModeTypeNames_
60  << exit(FatalError);
61 
62  return volumeModeType(0);
63 }
64 
65 
66 template<class Type>
68 (
69  const volumeModeType& vmtType
70 ) const
71 {
72  if (vmtType > volumeModeTypeNames_.size())
73  {
74  return "UNKNOWN";
75  }
76  else
77  {
78  return volumeModeTypeNames_[vmtType];
79  }
80 }
81 
82 
83 template<class Type>
85 {
86  fieldNames_.setSize(dict.toc().size());
87  injectionRate_.setSize(fieldNames_.size());
88 
89  applied_.setSize(fieldNames_.size(), false);
90 
91  label i = 0;
92  forAllConstIter(dictionary, dict, iter)
93  {
94  fieldNames_[i] = iter().keyword();
95  dict.lookup(iter().keyword()) >> injectionRate_[i];
96  i++;
97  }
98 
99  // Set volume normalisation
100  if (volumeMode_ == vmAbsolute)
101  {
102  VDash_ = V_;
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
109 template<class Type>
111 (
112  const word& name,
113  const word& modelType,
114  const dictionary& dict,
115  const fvMesh& mesh
116 )
117 :
118  cellSetOption(name, modelType, dict, mesh),
119  volumeMode_(vmAbsolute),
120  VDash_(1.0),
121  injectionRate_()
122 {
123  read(dict);
124 }
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
129 template<class Type>
131 (
132  fvMatrix<Type>& eqn,
133  const label fieldi
134 )
135 {
136  if (debug)
137  {
138  Info<< "SemiImplicitSource<" << pTraits<Type>::typeName
139  << ">::addSup for source " << name_ << endl;
140  }
141 
143 
145  (
146  IOobject
147  (
148  name_ + fieldNames_[fieldi] + "Su",
149  mesh_.time().timeName(),
150  mesh_,
151  IOobject::NO_READ,
152  IOobject::NO_WRITE
153  ),
154  mesh_,
156  (
157  "zero",
158  eqn.dimensions()/dimVolume,
159  Zero
160  ),
161  false
162  );
163 
164  UIndirectList<Type>(Su, cells_) = injectionRate_[fieldi].first()/VDash_;
165 
167  (
168  IOobject
169  (
170  name_ + fieldNames_[fieldi] + "Sp",
171  mesh_.time().timeName(),
172  mesh_,
173  IOobject::NO_READ,
174  IOobject::NO_WRITE
175  ),
176  mesh_,
178  (
179  "zero",
180  Su.dimensions()/psi.dimensions(),
181  0.0
182  ),
183  false
184  );
185 
186  UIndirectList<scalar>(Sp, cells_) = injectionRate_[fieldi].second()/VDash_;
187 
188  eqn += Su + fvm::SuSp(Sp, psi);
189 }
190 
191 
192 template<class Type>
194 (
195  const volScalarField& rho,
196  fvMatrix<Type>& eqn,
197  const label fieldi
198 )
199 {
200  if (debug)
201  {
202  Info<< "SemiImplicitSource<" << pTraits<Type>::typeName
203  << ">::addSup for source " << name_ << endl;
204  }
205 
206  return this->addSup(eqn, fieldi);
207 }
208 
209 
210 // ************************************************************************* //
volumeModeType wordToVolumeModeType(const word &vtName) const
Helper function to convert from a word to a volumeModeType.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
tmp< GeometricField< Type, fvPatchField, volMesh > > SuSp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:102
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:284
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Traits class for primitives.
Definition: pTraits.H:50
Generic GeometricField class.
Generic dimensioned Type class.
static const wordList volumeModeTypeNames_
Word list of volume mode type names.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
wordList toc() const
Return the table of contents.
Definition: dictionary.C:776
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
A class for handling words, derived from string.
Definition: word.H:59
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 volumeModeTypeToWord(const volumeModeType &vtType) const
Helper function to convert from a volumeModeType to a word.
static const zero Zero
Definition: zero.H:91
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
static const char nl
Definition: Ostream.H:262
volumeModeType
Enumeration for volume types.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
virtual void addSup(fvMatrix< Type > &eqn, const label fieldi)
Add explicit contribution to equation.
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 abtract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:69
const volScalarField & psi
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
SemiImplicitSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
void setFieldData(const dictionary &dict)
Set the local field data.
Calculate the matrix for implicit and explicit sources.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
zeroField Sp
Definition: alphaSuSp.H:2
const dimensionSet & dimensions() const
Definition: fvMatrix.H:289