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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2023 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 "patchInjection.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38 namespace surfaceFilmModels
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
51  const dictionary& dict
52 )
53 :
54  injectionModel(type(), film, dict),
55  deltaStable_(coeffDict_.getOrDefault<scalar>("deltaStable", 0))
56 {
58 
59  const label nPatches
60  (
61  pbm.size()
63  );
64 
66  if (coeffDict_.readIfPresent("patches", patchNames))
67  {
68  // Can also use pbm.indices(), but no warnings...
69  patchIDs_ = pbm.patchSet(patchNames).sortedToc();
70 
71  Info<< " applying to patches:" << nl;
72 
73  for (const label patchi : patchIDs_)
74  {
75  Info<< " " << pbm[patchi].name() << endl;
76  }
77  }
78  else
79  {
80  Info<< " applying to all patches" << endl;
81 
83  }
84 
86 
87  if (patchIDs_.empty())
88  {
90  << "No patches selected"
91  << exit(FatalError);
92  }
93 }
94 
95 
96 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
97 
99 (
100  scalarField& availableMass,
101  scalarField& massToInject,
102  scalarField& diameterToInject
103 )
104 {
105  // Do not correct if no patches selected
106  if (patchIDs_.empty()) return;
107 
108  const scalarField& delta = film().delta();
109  const scalarField& rho = film().rho();
110  const scalarField& magSf = film().magSf();
111 
113 
114  forAll(patchIDs_, pidi)
115  {
116  label patchi = patchIDs_[pidi];
117  const polyPatch& pp = pbm[patchi];
118  const labelList& faceCells = pp.faceCells();
119 
120  // Accumulate the total mass removed from patch
121  scalar dMassPatch = 0;
122 
123  for (const label celli : faceCells)
124  {
125  scalar ddelta = max(0.0, delta[celli] - deltaStable_);
126  scalar dMass = ddelta*rho[celli]*magSf[celli];
127  massToInject[celli] += dMass;
128  availableMass[celli] -= dMass;
129  dMassPatch += dMass;
130  }
131 
132  patchInjectedMasses_[pidi] += dMassPatch;
133  addToInjectedMass(dMassPatch);
134  }
135 
137 
138  if (writeTime())
139  {
140  scalarField patchInjectedMasses0
141  (
142  getModelProperty<scalarField>
143  (
144  "patchInjectedMasses",
146  )
147  );
148 
151  patchInjectedMasses0 += patchInjectedMassTotals;
152 
153  setModelProperty<scalarField>
154  (
155  "patchInjectedMasses",
156  patchInjectedMasses0
157  );
158 
160  }
161 }
162 
163 
165 {
166  // Do not correct if no patches selected
167  if (!patchIDs_.size()) return;
168 
169  scalarField patchInjectedMasses
170  (
171  getModelProperty<scalarField>
172  (
173  "patchInjectedMasses",
175  )
176  );
177 
180 
181  forAll(patchIDs_, pidi)
182  {
183  label patchi = patchIDs_[pidi];
184  patchMasses[patchi] +=
185  patchInjectedMasses[pidi] + patchInjectedMassTotals[pidi];
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 } // End namespace surfaceFilmModels
193 } // End namespace regionModels
194 } // End namespace Foam
195 
196 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
scalar delta
const polyBoundaryMesh & pbm
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
patchInjection(const patchInjection &)=delete
No copy construct.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
scalarField patchInjectedMasses_
Injected mass for each patch at which the film is removed.
virtual void patchInjectedMassTotals(scalarField &patchMasses) const
Accumulate the total mass injected for the patches into the.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Remove and inject the mass in the film as it passes over the selected patches.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:92
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
virtual bool writeTime() const
Flag to indicate when to write a property.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
labelList patchIDs_
List of patch IDs at which the film is removed.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList patchNames(nPatches)
virtual const volScalarField & delta() const =0
Return the film thickness [m].
scalar deltaStable_
Stable film thickness - mass only removed if thickness exceeds.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void addToInjectedMass(const scalar dMass)
Add to injected mass.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Base class for film injection models, handling mass transfer from the film.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
messageStream Info
Information stream (stdout output on master, null elsewhere)
static void listCombineGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
defineTypeNameAndDebug(kinematicSingleLayer, 0)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127