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 -------------------------------------------------------------------------------
11 License
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 
33 namespace Foam
34 {
35 namespace regionModels
36 {
37 namespace surfaceFilmModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 injectionModelList::injectionModelList(surfaceFilmRegionModel& film)
43 :
45  filmSubModelBase(film)
46 {}
47 
48 
49 injectionModelList::injectionModelList
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 // ************************************************************************* //
setSize
points setSize(newPointi)
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::regionModels::surfaceFilmModels::injectionModel::correct
void correct()
Correct.
Definition: injectionModel.C:81
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
injectionModelList.H
Foam::subModelBase::getBaseProperty
Type getBaseProperty(const word &entryName, const Type &defaultValue=Type(Zero)) const
Retrieve generic property from the base model.
Definition: subModelBaseTemplates.C:32
Foam::regionModels::surfaceFilmModels::filmSubModelBase::writeTime
virtual bool writeTime() const
Flag to indicate when to write a property.
Definition: filmSubModelBase.C:98
Foam::regionModels::surfaceFilmModels::injectionModelList::info
virtual void info(Ostream &os)
Provide some info.
Definition: injectionModelList.C:126
Foam::regionModels::surfaceFilmModels::filmSubModelBase::film
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
Definition: filmSubModelBaseI.H:39
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::HashSet< word, Hash< word > >
Foam::regionModels::surfaceFilmModels::injectionModel::New
static autoPtr< injectionModel > New(surfaceFilmRegionModel &film, const dictionary &dict, const word &mdoelType)
Return a reference to the selected injection model.
Definition: injectionModelNew.C:43
Foam::globalMeshData::processorPatches
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
Definition: globalMeshData.H:381
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::regionModels::surfaceFilmModels::injectionModel::patchInjectedMassTotals
virtual void patchInjectedMassTotals(scalarField &patchMasses) const
Accumulate the total mass injected for the patches into the.
Definition: injectionModel.H:153
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
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:123
Foam::regionModels::surfaceFilmModels::injectionModel
Base class for film injection models, handling mass transfer from the film.
Definition: injectionModel.H:58
os
OBJstream os(runTime.globalPath()/outputName)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:339
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
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:290
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::regionModels::surfaceFilmModels::injectionModelList::correct
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Definition: injectionModelList.C:98
Foam::List< word >
Foam::subModelBase::setBaseProperty
void setBaseProperty(const word &entryName, const Type &value)
Add generic property to the base model.
Definition: subModelBaseTemplates.C:66
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::regionModels::surfaceFilmModels::injectionModel::injectedMassTotal
virtual scalar injectedMassTotal() const
Return the total mass injected.
Definition: injectionModel.C:93
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::plusEqOp
Definition: ops.H:72
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1295
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::regionModels::regionModel::intCoupledPatchIDs
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:175
Foam::regionModels::surfaceFilmModels::filmSubModelBase
Base class for surface film sub-models.
Definition: filmSubModelBase.H:56
Foam::regionModels::surfaceFilmModels::injectionModelList::~injectionModelList
virtual ~injectionModelList()
Destructor.
Definition: injectionModelList.C:91
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62