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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "injectionModelList.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace regionModels
36{
37namespace surfaceFilmModels
38{
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
43:
46{}
47
48
50(
52 const dictionary& dict
53)
54:
57 (
58 "injectionModelList",
59 film,
60 dict,
61 "injectionModelList",
62 "injectionModelList"
63 ),
64 massInjected_(film.intCoupledPatchIDs().size(), 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 * * * * * * * * * * * * * * * //
90
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 forAll(*this, i)
106 {
107 injectionModel& im = operator[](i);
108 im.correct(availableMass, massToInject, diameterToInject);
109 }
110
111 // Push values to boundaries ready for transfer to the primary region
112 massToInject.correctBoundaryConditions();
113 diameterToInject.correctBoundaryConditions();
114
115 const labelList& patchIDs = film().intCoupledPatchIDs();
116
117 forAll(patchIDs, i)
118 {
119 label patchi = patchIDs[i];
120 massInjected_[i] =
121 massInjected_[i] + sum(massToInject.boundaryField()[patchi]);
122 }
123}
124
125
127{
128 const polyBoundaryMesh& pbm = film().regionMesh().boundaryMesh();
129
130 scalar injectedMass = 0;
131 scalarField patchInjectedMasses
132 (
133 pbm.size() - film().regionMesh().globalData().processorPatches().size(),
134 0
135 );
136
137 forAll(*this, i)
138 {
139 const injectionModel& im = operator[](i);
140 injectedMass += im.injectedMassTotal();
141 im.patchInjectedMassTotals(patchInjectedMasses);
142 }
143
144 os << indent << "injected mass = " << injectedMass << nl;
145
146 forAll(patchInjectedMasses, patchi)
147 {
148 if (mag(patchInjectedMasses[patchi]) > VSMALL)
149 {
150 os << indent << indent << "from patch " << pbm[patchi].name()
151 << " = " << patchInjectedMasses[patchi] << nl;
152 }
153 }
154
155 scalarField mass0(massInjected_.size(), Zero);
156 this->getBaseProperty("massInjected", mass0);
157
158 scalarField mass(massInjected_);
160 mass += mass0;
161
162 const labelList& patchIDs = film().intCoupledPatchIDs();
163
164 forAll(patchIDs, i)
165 {
166 label patchi = patchIDs[i];
167 Info<< indent << " - patch: " << pbm[patchi].name() << ": "
168 << mass[i] << endl;
169 }
170
171 if (film().time().writeTime())
172 {
173 setBaseProperty("massInjected", mass);
174 massInjected_ = 0.0;
175 }
176}
177
178
179// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180
181} // End namespace surfaceFilmModels
182} // End namespace regionModels
183} // End namespace Foam
184
185// ************************************************************************* //
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static void listCombineGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const injectionModel * set(const label i) const
Definition: PtrList.H:138
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:151
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const T & operator[](const label i) const
Return const reference to the element.
Definition: UPtrListI.H:234
label size() const noexcept
The number of elements 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:126
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
InfoProxy< ensightCells > info() const
Return info proxy.
Definition: ensightCells.H:254
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:175
Base class for surface film sub-models.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
virtual bool writeTime() const
Flag to indicate when to write a property.
Base class for film injection models, handling mass transfer from the film.
virtual scalar injectedMassTotal() const
Return the total mass injected.
virtual void patchInjectedMassTotals(scalarField &patchMasses) const
Accumulate the total mass injected for the patches into the.
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.
Definition: subModelBase.C:113
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:68
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333