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 | Copyright (C) 2011-2015 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 "injectionModelList.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace regionModels
33 {
34 namespace surfaceFilmModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 :
42  filmSubModelBase(owner)
43 {}
44 
45 
47 (
48  surfaceFilmModel& owner,
49  const dictionary& dict
50 )
51 :
54  (
55  "injectionModelList",
56  owner,
57  dict,
58  "injectionModelList",
59  "injectionModelList"
60  ),
61  massInjected_(owner.intCoupledPatchIDs().size(), 0.0)
62 {
63  const wordList activeModels(dict.lookup("injectionModels"));
64 
65  wordHashSet models;
66  forAll(activeModels, i)
67  {
68  models.insert(activeModels[i]);
69  }
70 
71  Info<< " Selecting film injection models" << endl;
72  if (models.size() > 0)
73  {
74  this->setSize(models.size());
75 
76  label i = 0;
77  forAllConstIter(wordHashSet, models, iter)
78  {
79  const word& model = iter.key();
80  set(i, injectionModel::New(owner, dict, model));
81  i++;
82  }
83  }
84  else
85  {
86  Info<< " none" << endl;
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 (
101  scalarField& availableMass,
102  volScalarField& massToInject,
103  volScalarField& diameterToInject
104 )
105 {
106  // Correct models that accumulate mass and diameter transfers
107  forAll(*this, i)
108  {
109  injectionModel& im = operator[](i);
110  im.correct(availableMass, massToInject, diameterToInject);
111  }
112 
113  // Push values to boundaries ready for transfer to the primary region
114  massToInject.correctBoundaryConditions();
115  diameterToInject.correctBoundaryConditions();
116 
117  const labelList& patchIDs = owner().intCoupledPatchIDs();
118 
119  forAll(patchIDs, i)
120  {
121  label patchi = patchIDs[i];
122  massInjected_[i] =
123  massInjected_[i] + sum(massToInject.boundaryField()[patchi]);
124  }
125 }
126 
127 
129 {
130  const polyBoundaryMesh& pbm = owner().regionMesh().boundaryMesh();
131 
132  scalar injectedMass = 0;
133  scalarField patchInjectedMasses(pbm.size(), 0);
134 
135  forAll(*this, i)
136  {
137  const injectionModel& im = operator[](i);
138  injectedMass += im.injectedMassTotal();
139  im.patchInjectedMassTotals(patchInjectedMasses);
140  }
141 
142  os << indent << "injected mass = " << injectedMass << nl;
143 
144  forAll(pbm, patchi)
145  {
146  if (mag(patchInjectedMasses[patchi]) > VSMALL)
147  {
148  os << indent << indent << "from patch " << pbm[patchi].name()
149  << " = " << patchInjectedMasses[patchi] << nl;
150  }
151  }
152 
153  scalarField mass0(massInjected_.size(), 0.0);
154  this->getBaseProperty("massInjected", mass0);
155 
158  mass += mass0;
159 
160  const labelList& patchIDs = owner().intCoupledPatchIDs();
161 
162  forAll(patchIDs, i)
163  {
164  label patchi = patchIDs[i];
165  Info<< indent << " - patch: " << pbm[patchi].name() << ": "
166  << mass[i] << endl;
167  }
168 
169  if (owner().time().outputTime())
170  {
171  setBaseProperty("massInjected", mass);
172  massInjected_ = 0.0;
173  }
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace surfaceFilmModels
180 } // End namespace regionModels
181 } // End namespace Foam
182 
183 // ************************************************************************* //
setSize
points setSize(newPointi)
injectionModelList.H
Foam::regionModels::surfaceFilmModels::injectionModel::correct
void correct()
Correct.
Definition: injectionModel.C:79
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::regionModels::surfaceFilmModels::injectionModel::New
static autoPtr< injectionModel > New(surfaceFilmModel &owner, const dictionary &dict, const word &mdoelType)
Return a reference to the selected injection model.
Definition: injectionModelNew.C:40
Foam::regionModels::surfaceFilmModels::injectionModelList::info
virtual void info(Ostream &os)
Provide some info.
Definition: injectionModelList.C:128
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::regionModels::surfaceFilmModels::filmSubModelBase::owner
const surfaceFilmModel & owner() const
Return const access to the owner surface film model.
Definition: filmSubModelBaseI.H:37
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::regionModels::surfaceFilmModels::filmSubModelBase::outputTime
virtual bool outputTime() const
Flag to indicate when to write a property.
Definition: filmSubModelBase.C:96
Foam::regionModels::surfaceFilmModels::injectionModel::patchInjectedMassTotals
virtual void patchInjectedMassTotals(scalarField &patchMasses) const
Accumulate the total mass injected for the patches into the.
Definition: injectionModel.H:151
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::PtrList< injectionModel >::operator[]
const T & operator[](const label) const
Return element const reference.
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::regionModels::surfaceFilmModels::injectionModelList::massInjected_
scalarField massInjected_
List of mass injected per patch.
Definition: injectionModelList.H:64
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Foam::regionModels::surfaceFilmModels::injectionModel
Base class for film injection models, handling mass transfer from the film.
Definition: injectionModel.H:56
Foam::plusEqOp
Definition: ops.H:71
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:282
Foam::subModelBase::getBaseProperty
Type getBaseProperty(const word &entryName, const Type &defaultValue=pTraits< Type >::zero) const
Retrieve generic property from the base model.
Definition: subModelBaseTemplates.C:30
Foam::regionModels::surfaceFilmModels::injectionModelList::correct
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Definition: injectionModelList.C:100
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::subModelBase::setBaseProperty
void setBaseProperty(const word &entryName, const Type &value)
Add generic property to the base model.
Definition: subModelBaseTemplates.C:64
Foam::regionModels::surfaceFilmModels::injectionModel::injectedMassTotal
virtual scalar injectedMassTotal() const
Return the total mass injected.
Definition: injectionModel.C:91
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::regionModels::regionModel::intCoupledPatchIDs
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:179
Foam::regionModels::surfaceFilmModels::filmSubModelBase
Definition: filmSubModelBase.H:54
Foam::regionModels::surfaceFilmModels::injectionModelList::~injectionModelList
virtual ~injectionModelList()
Destructor.
Definition: injectionModelList.C:93
Foam::regionModels::surfaceFilmModels::injectionModelList::injectionModelList
injectionModelList(const injectionModelList &)
Disallow default bitwise copy construct.
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61