patchInjection.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) 2015-2016 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 "patchInjection.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace regionModels
34 {
35 namespace surfaceFilmModels
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 defineTypeNameAndDebug(patchInjection, 0);
41 addToRunTimeSelectionTable(injectionModel, patchInjection, dictionary);
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 patchInjection::patchInjection
46 (
47  surfaceFilmModel& owner,
48  const dictionary& dict
49 )
50 :
51  injectionModel(type(), owner, dict),
52  deltaStable_(coeffDict_.lookupOrDefault<scalar>("deltaStable", 0.0))
53 {
54  const polyBoundaryMesh& pbm = owner.regionMesh().boundaryMesh();
55  patchIDs_.setSize(pbm.size());
56 
57  if (coeffDict_.found("patches"))
58  {
59  wordReList patchNames(coeffDict_.lookup("patches"));
60  const labelHashSet patchSet = pbm.patchSet(patchNames);
61 
62  Info<< " applying to patches:" << nl;
63 
64  label pidi = 0;
65  forAllConstIter(labelHashSet, patchSet, iter)
66  {
67  label patchi = iter.key();
68  patchIDs_[pidi++] = patchi;
69  Info<< " " << pbm[patchi].name() << endl;
70  }
71  patchIDs_.setSize(pidi);
72  patchInjectedMasses_.setSize(pidi, 0);
73  }
74  else
75  {
76  Info<< " applying to all patches" << endl;
77 
78  forAll(patchIDs_, patchi)
79  {
80  patchIDs_[patchi] = patchi;
81  }
82 
83  patchInjectedMasses_.setSize(pbm.size(), 0);
84  }
85 
86  if (!patchIDs_.size())
87  {
89  << "No patches selected"
90  << exit(FatalError);
91  }
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
98 {}
99 
100 
101 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
102 
104 (
105  scalarField& availableMass,
106  scalarField& massToInject,
107  scalarField& diameterToInject
108 )
109 {
110  // Do not correct if no patches selected
111  if (!patchIDs_.size()) return;
112 
113  const scalarField& delta = owner().delta();
114  const scalarField& rho = owner().rho();
115  const scalarField& magSf = owner().magSf();
116 
117  const polyBoundaryMesh& pbm = owner().regionMesh().boundaryMesh();
118 
119  forAll(patchIDs_, pidi)
120  {
121  label patchi = patchIDs_[pidi];
122  const polyPatch& pp = pbm[patchi];
123  const labelList& faceCells = pp.faceCells();
124 
125  // Accumulate the total mass removed from patch
126  scalar dMassPatch = 0;
127 
128  forAll(faceCells, fci)
129  {
130  label celli = faceCells[fci];
131 
132  scalar ddelta = max(0.0, delta[celli] - deltaStable_);
133  scalar dMass = ddelta*rho[celli]*magSf[celli];
134  massToInject[celli] += dMass;
135  availableMass[celli] -= dMass;
136  dMassPatch += dMass;
137  }
138 
139  patchInjectedMasses_[pidi] += dMassPatch;
140  addToInjectedMass(dMassPatch);
141  }
142 
144 
145  if (writeTime())
146  {
147  scalarField patchInjectedMasses0
148  (
149  getModelProperty<scalarField>
150  (
151  "patchInjectedMasses",
153  )
154  );
155 
157  Pstream::listCombineGather(patchInjectedMassTotals, plusEqOp<scalar>());
158  patchInjectedMasses0 += patchInjectedMassTotals;
159 
160  setModelProperty<scalarField>
161  (
162  "patchInjectedMasses",
163  patchInjectedMasses0
164  );
165 
167  }
168 }
169 
170 
172 {
173  // Do not correct if no patches selected
174  if (!patchIDs_.size()) return;
175 
176  scalarField patchInjectedMasses
177  (
178  getModelProperty<scalarField>
179  (
180  "patchInjectedMasses",
182  )
183  );
184 
186  Pstream::listCombineGather(patchInjectedMassTotals, plusEqOp<scalar>());
187 
188  forAll(patchIDs_, pidi)
189  {
190  label patchi = patchIDs_[pidi];
191  patchMasses[patchi] +=
192  patchInjectedMasses[pidi] + patchInjectedMassTotals[pidi];
193  }
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace surfaceFilmModels
200 } // End namespace regionModels
201 } // End namespace Foam
202 
203 // ************************************************************************* //
scalar delta
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
scalarField patchInjectedMasses_
Injected mass for each patch at which the film is removed.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual void patchInjectedMassTotals(scalarField &patchMasses) const
Accumulate the total mass injected for the patches into the.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Macros for easy insertion into run-time selection tables.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
labelList patchIDs_
List of patch IDs at which the film is removed.
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList patchNames(nPatches)
scalar deltaStable_
Stable film thickness - mass only removed if thickness execeeds.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual const volScalarField & delta() const =0
Return the film thickness [m].
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
void addToInjectedMass(const scalar dMass)
Add to injected mass.
static const char nl
Definition: Ostream.H:262
Base class for film injection models, handling mass transfer from the film.
label patchi
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
const surfaceFilmModel & owner() const
Return const access to the owner surface film model.
messageStream Info
virtual bool writeTime() const
Flag to indicate when to write a property.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
defineTypeNameAndDebug(kinematicSingleLayer, 0)
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].