Loading...
Searching...
No Matches
injectionModelList.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) 2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "injectionModelList.H"
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace regionModels
36{
38{
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
42injectionModelList::injectionModelList(liquidFilmBase& film)
43:
46{}
47
48
49injectionModelList::injectionModelList
50(
51 liquidFilmBase& film,
52 const dictionary& dict
53)
54:
57 (
58 "injectionModelList",
59 film,
60 dict,
61 "injectionModelList",
62 "injectionModelList"
63 ),
64 massInjected_(Zero)
65{
66 const wordList activeModels(dict.lookup("injectionModels"));
67
68 wordHashSet models(activeModels);
69
70 Info<< " Selecting film injection models" << endl;
71 if (models.size())
72 {
73 this->setSize(models.size());
74
75 label i = 0;
76 for (const word& model : models)
77 {
78 set(i, injectionModel::New(film, dict, model));
79 i++;
80 }
81 }
82 else
83 {
84 Info<< " none" << endl;
85 }
86}
87
88
89// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92{}
93
94
95// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96
98(
99 scalarField& availableMass,
100 volScalarField& massToInject,
101 volScalarField& diameterToInject
102)
103{
104 // Correct models that accumulate mass and diameter transfers
105 PtrList<injectionModel>& models = *this;
107
108 if (models.empty() || patchIds.empty())
109 {
110 // Nothing to do
111 return;
112 }
113
114 const label patchi = patchIds[0];
115
116 if
117 (
118 patchIds.size() == 1
119 &&
120 (
121 film().regionMesh().nFaces()
122 == massToInject.boundaryField()[patchi].size()
123 )
124 )
125 {
126 // All faces of a single patch.
127 //- Can write directly into boundary fiels
128
129 for (injectionModel& im : models)
130 {
131 im.correct
132 (
133 availableMass,
134 massToInject.boundaryFieldRef()[patchi],
135 diameterToInject.boundaryFieldRef()[patchi]
136 );
137 }
138
139 massInjected_ += gSum(massToInject.boundaryField()[patchi]);
140 }
141 else
142 {
143 // Multiple patches, or a subset of faces
144 // - probably needs intermediate variables
145
146 // The polyPatch/local-face for each of the faceLabels
147 const List<labelPair>& patchFaces
149
150 if (patchFaces.size() != availableMass.size())
151 {
153 << "film has " << patchFaces.size()
154 << " faces, but availableMass has " << availableMass.size()
155 << " entries"
156 << abort(FatalError);
157 }
158
159 scalarField massToInjectTmp(availableMass.size(), Zero);
160 scalarField diameterToInjectTmp(availableMass.size(), Zero);
161
162 for (injectionModel& im : models)
163 {
164 im.correct
165 (
166 availableMass,
167 massToInjectTmp,
168 diameterToInjectTmp
169 );
170 }
171
172 massInjected_ += gSum(massToInjectTmp);
173
174 // Transcribe to boundary patches
175 forAll(patchFaces, i)
176 {
177 const labelPair& patchAndFace = patchFaces[i];
178
179 if (patchAndFace.first() >= 0)
180 {
181 massToInject.boundaryFieldRef()[patchAndFace]
182 = massToInjectTmp[i];
183
184 diameterToInject.boundaryFieldRef()[patchAndFace]
185 = diameterToInjectTmp[i];
186 }
188 }
189}
190
191
193{
195
196 scalar injectedMass = 0;
197 scalar patchInjectedMasses = 0;
198
199 for (const injectionModel& im : *this)
200 {
201 injectedMass += im.injectedMassTotal();
202 im.patchInjectedMassTotals(patchInjectedMasses);
203 }
204
205 scalar mass0(Zero);
206 this->getBaseProperty("massInjected", mass0);
207
208 scalar mass(massInjected_);
209 mass += mass0;
210
211 os << indent << "injected mass = " << injectedMass << nl;
212
213 if (mag(patchInjectedMasses) > VSMALL)
214 {
215 os << indent << indent << "from patch ";
216
217 for (const label patchi : film().primaryPatchIDs())
218 {
219 os << ' ' << pbm[patchi].name();
220 }
221 os << " = " << patchInjectedMasses << nl;
222 }
223
224 {
225 Info<< indent << " - patch:";
226 for (const label patchi : film().primaryPatchIDs())
227 {
228 os << ' ' << pbm[patchi].name();
229 }
230 os << " " << mass << endl;
231 }
232
233 if (film().primaryMesh().time().writeTime())
234 {
235 setBaseProperty("massInjected", mass);
236 massInjected_ = 0.0;
237 }
238}
239
240
241// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242
243} // End namespace areaSurfaceFilmModels
244} // End namespace regionModels
245} // End namespace Foam
246
247// ************************************************************************* //
const polyBoundaryMesh & pbm
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const injectionModel * set(const label i) const
Definition PtrList.H:171
void setSize(const label n)
Definition PtrList.H:357
constexpr PtrList() noexcept
Definition PtrListI.H:29
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
bool empty() const noexcept
True if the list is empty (ie, size() is zero).
Definition UPtrListI.H:99
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels().
Definition faMeshI.H:138
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const liquidFilmBase & film() const
Return const access to the film surface film model.
virtual bool writeTime() const
Flag to indicate when to write a property.
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Base class for film injection models, handling mass transfer from the film.
virtual scalar injectedMassTotal() const
Return the total mass injected.
virtual void patchInjectedMassTotals(scalar &patchMasses) const
Accumulate the total mass injected for the patches into the scalarField provided.
static autoPtr< injectionModel > New(liquidFilmBase &film, const dictionary &dict, const word &mdoelType)
Return a reference to the selected injection model.
const faMesh & regionMesh() const
Return the region mesh database.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
const labelList & primaryPatchIDs() const
List of patch IDs on the primary region coupled to this region.
Type getBaseProperty(const word &entryName, const Type &defaultValue=Type(Zero)) const
Retrieve generic property from the base model.
const dictionary & dict() const
Return const access to the cloud dictionary.
void setBaseProperty(const word &entryName, const Type &value)
Add generic property to the base model.
A class for handling words, derived from Foam::string.
Definition word.H:66
labelList patchIds
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
Type gSum(const FieldField< Field, Type > &f)
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299