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 -------------------------------------------------------------------------------
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 "injectionModelList.H"
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace regionModels
36 {
37 namespace areaSurfaceFilmModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 injectionModelList::injectionModelList(liquidFilmBase& film)
43 :
45  filmSubModelBase(film)
46 {}
47 
48 
49 injectionModelList::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 * * * * * * * * * * * * * * * //
90 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 (
99  scalarField& availableMass,
100  volScalarField& massToInject,
101  volScalarField& diameterToInject
102 )
103 {
104  const label patchi = film().patchID();
105 
106  // Correct models that accumulate mass and diameter transfers
107  forAll(*this, i)
108  {
109  injectionModel& im = operator[](i);
110  im.correct
111  (
112  availableMass,
113  massToInject.boundaryFieldRef()[patchi],
114  diameterToInject.boundaryFieldRef()[patchi]
115  );
116  }
117 
118  massInjected_ += gSum(massToInject.boundaryField()[patchi]);
119 }
120 
121 
123 {
124  const polyBoundaryMesh& pbm = film().primaryMesh().boundaryMesh();
125 
126  scalar injectedMass = 0;
127  scalar patchInjectedMasses = 0;
128 
129  forAll(*this, i)
130  {
131  const injectionModel& im = operator[](i);
132  injectedMass += im.injectedMassTotal();
133  im.patchInjectedMassTotals(patchInjectedMasses);
134  }
135 
136  os << indent << "injected mass = " << injectedMass << nl;
137 
138  const label patchi = film().patchID();
139 
140  if (mag(patchInjectedMasses) > VSMALL)
141  {
142  os << indent << indent << "from patch " << pbm[patchi].name()
143  << " = " << patchInjectedMasses << nl;
144  }
145 
146  scalar mass0(Zero);
147  this->getBaseProperty("massInjected", mass0);
148 
149  scalar mass(massInjected_);
150  mass += mass0;
151 
152  Info<< indent << " - patch: " << pbm[patchi].name() << " "
153  << mass << endl;
154 
155  if (film().primaryMesh().time().writeTime())
156  {
157  setBaseProperty("massInjected", mass);
158  massInjected_ = 0.0;
159  }
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace areaSurfaceFilmModels
166 } // End namespace regionModels
167 } // End namespace Foam
168 
169 // ************************************************************************* //
injectionModelList.H
volFields.H
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase
Definition: liquidFilmBase.H:56
setSize
points setSize(newPointi)
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Definition: BitOps.C:30
Foam::regionModels::areaSurfaceFilmModels::injectionModel::correct
void correct()
Definition: injectionModel.C:74
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::regionModels::areaSurfaceFilmModels::injectionModelList::~injectionModelList
virtual ~injectionModelList()
Definition: injectionModelList.C:84
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:59
Foam::regionModels::regionFaModel::primaryMesh
const fvMesh & primaryMesh() const
Definition: regionFaModelI.H:26
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::subModelBase::getBaseProperty
Type getBaseProperty(const word &entryName, const Type &defaultValue=Type(Zero)) const
Definition: subModelBaseTemplates.C:25
Foam::regionModels::areaSurfaceFilmModels::injectionModelList::info
virtual void info(Ostream &os)
Definition: injectionModelList.C:115
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::PtrList< injectionModel >::set
const T * set(const label i) const
Definition: PtrList.H:134
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:587
Foam::regionModels::areaSurfaceFilmModels::injectionModel::patchInjectedMassTotals
virtual void patchInjectedMassTotals(scalar &patchMasses) const
Definition: injectionModel.H:149
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:73
Foam::regionModels::areaSurfaceFilmModels::filmSubModelBase::writeTime
virtual bool writeTime() const
Definition: filmSubModelBase.C:91
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::regionModels::areaSurfaceFilmModels::injectionModelList::correct
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Definition: injectionModelList.C:91
Foam::regionModels::areaSurfaceFilmModels::injectionModel
Base class for film injection models, handling mass transfer from the film.
Definition: injectionModel.H:54
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::subModelBase::dict
const dictionary & dict() const
Definition: subModelBase.C:106
Foam::Info
messageStream Info
Foam::regionModels::areaSurfaceFilmModels::filmSubModelBase::film
const liquidFilmBase & film() const
Definition: filmSubModelBaseI.H:32
Foam::PtrList< injectionModel >::setSize
void setSize(const label newLen)
Definition: PtrList.H:147
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:379
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
os
OBJstream os(runTime.globalPath()/outputName)
Foam
Definition: atmBoundaryLayer.C:26
Foam::indent
Ostream & indent(Ostream &os)
Definition: Ostream.H:343
Foam::regionModels::areaSurfaceFilmModels::filmSubModelBase
Definition: filmSubModelBase.H:52
Foam::IOobject::name
const word & name() const noexcept
Definition: IOobjectI.H:58
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Definition: GeometricField.C:776
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::subModelBase::setBaseProperty
void setBaseProperty(const word &entryName, const Type &value)
Definition: subModelBaseTemplates.C:59
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::regionModels::areaSurfaceFilmModels::injectionModel::injectedMassTotal
virtual scalar injectedMassTotal() const
Definition: injectionModel.C:86
Foam::regionModels::regionFaModel::patchID
label patchID() const
Definition: regionFaModelI.H:131
Foam::regionModels::areaSurfaceFilmModels::injectionModel::New
static autoPtr< injectionModel > New(liquidFilmBase &film, const dictionary &dict, const word &mdoelType)
Definition: injectionModelNew.C:35
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Definition: GeometricFieldI.H:55