transferModelList.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) 2017 OpenFOAM Foundation
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 "transferModelList.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace surfaceFilmModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 transferModelList::transferModelList(surfaceFilmRegionModel& film)
42 :
44  filmSubModelBase(film)
45 {}
46 
47 
48 transferModelList::transferModelList
49 (
51  const dictionary& dict
52 )
53 :
56  (
57  "transferModelList",
58  film,
59  dict,
60  "transferModelList",
61  "transferModelList"
62  ),
63  massTransferred_(film.intCoupledPatchIDs().size(), Zero)
64 {
65  const wordList activeModels
66  (
67  dict.lookupOrDefault("transferModels", wordList())
68  );
69 
70  wordHashSet models(activeModels);
71 
72  Info<< " Selecting film transfer models" << endl;
73  if (models.size() > 0)
74  {
75  this->setSize(models.size());
76 
77  label i = 0;
78  for (const word& model : models)
79  {
80  set(i, transferModel::New(film, 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& massToTransfer
103 )
104 {
105  // Correct models that accumulate mass and diameter transfers
106  forAll(*this, i)
107  {
108  operator[](i).correct(availableMass, massToTransfer);
109  }
110 
111  // Push values to boundaries ready for transfer to the primary region
112  massToTransfer.correctBoundaryConditions();
113 
114  const labelList& patchIDs = film().intCoupledPatchIDs();
115 
116  forAll(patchIDs, i)
117  {
118  label patchi = patchIDs[i];
119  massTransferred_[i] =
120  massTransferred_[i] + sum(massToTransfer.boundaryField()[patchi]);
121  }
122 }
123 
124 
126 (
127  scalarField& availableMass,
128  volScalarField& massToTransfer,
129  volScalarField& energyToTransfer
130 )
131 {
132  // Correct models that accumulate mass and diameter transfers
133  forAll(*this, i)
134  {
135  operator[](i).correct(availableMass, massToTransfer, energyToTransfer);
136  }
137 
138  // Push values to boundaries ready for transfer to the primary region
139  massToTransfer.correctBoundaryConditions();
140  energyToTransfer.correctBoundaryConditions();
141 
142  const labelList& patchIDs = film().intCoupledPatchIDs();
143 
144  forAll(patchIDs, i)
145  {
146  label patchi = patchIDs[i];
147  massTransferred_[i] =
148  massTransferred_[i] + sum(massToTransfer.boundaryField()[patchi]);
149  }
150 }
151 
152 
154 {
155  const polyBoundaryMesh& pbm = film().regionMesh().boundaryMesh();
156 
157  scalar transferredMass = 0;
158  scalarField patchTransferredMasses
159  (
160  pbm.size() - film().regionMesh().globalData().processorPatches().size(),
161  0
162  );
163 
164  forAll(*this, i)
165  {
166  const transferModel& im = operator[](i);
167  transferredMass += im.transferredMassTotal();
168  im.patchTransferredMassTotals(patchTransferredMasses);
169  }
170 
171  os << indent << "transferred mass = " << transferredMass << nl;
172 
173  forAll(patchTransferredMasses, patchi)
174  {
175  if (mag(patchTransferredMasses[patchi]) > VSMALL)
176  {
177  os << indent << indent << "from patch " << pbm[patchi].name()
178  << " = " << patchTransferredMasses[patchi] << nl;
179  }
180  }
181 
182  scalarField mass0(massTransferred_.size(), Zero);
183  this->getBaseProperty("massTransferred", mass0);
184 
185  scalarField mass(massTransferred_);
187  mass += mass0;
188 
189  const labelList& patchIDs = film().intCoupledPatchIDs();
190 
191  forAll(patchIDs, i)
192  {
193  label patchi = patchIDs[i];
194  Info<< indent << " - patch: " << pbm[patchi].name() << ": "
195  << mass[i] << endl;
196  }
197 
198  if (film().time().writeTime())
199  {
200  setBaseProperty("massTransferred", mass);
201  massTransferred_ = 0.0;
202  }
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 } // End namespace surfaceFilmModels
209 } // End namespace regionModels
210 } // End namespace Foam
211 
212 // ************************************************************************* //
setSize
points setSize(newPointi)
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:46
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::regionModels::surfaceFilmModels::transferModelList::correct
virtual void correct(scalarField &availableMass, volScalarField &massToTransfer)
Correct kinematic transfers.
Definition: transferModelList.C:100
Foam::regionModels::surfaceFilmModels::transferModelList::info
virtual void info(Ostream &os)
Provide some info.
Definition: transferModelList.C:153
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::globalMeshData::processorPatches
const labelList & processorPatches() const
Return list of processor patch labels.
Definition: globalMeshData.H:404
transferModelList.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::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:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::HashSet< word >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::regionModels::surfaceFilmModels::transferModel::transferredMassTotal
virtual scalar transferredMassTotal() const
Return the total mass transferred.
Definition: transferModel.C:105
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:62
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:63
Foam::regionModels::surfaceFilmModels::transferModel::patchTransferredMassTotals
virtual void patchTransferredMassTotals(scalarField &patchMasses) const
Accumulate the total mass transferred for the patches into the.
Definition: transferModel.H:160
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
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:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:307
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:909
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::regionModels::surfaceFilmModels::transferModel
Base class for film transfer models, handling mass transfer between the film and the continuous phase...
Definition: transferModel.H:58
Foam::nl
constexpr char nl
Definition: Ostream.H:372
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::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::regionModels::surfaceFilmModels::transferModelList::~transferModelList
virtual ~transferModelList()
Destructor.
Definition: transferModelList.C:93
Foam::regionModels::surfaceFilmModels::transferModel::New
static autoPtr< transferModel > New(surfaceFilmRegionModel &film, const dictionary &dict, const word &modelType)
Return a reference to the selected injection model.
Definition: transferModelNew.C:43
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:1241
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:182
Foam::regionModels::surfaceFilmModels::filmSubModelBase
Base class for surface film sub-models.
Definition: filmSubModelBase.H:56
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62