Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
34namespace Foam
35{
36namespace regionModels
38namespace 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{
57 const polyBoundaryMesh& pbm = film.regionMesh().boundaryMesh();
58
59 const label nPatches
60 (
61 pbm.size()
62 - film.regionMesh().globalData().processorPatches().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"
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
119 // Accumulate the total mass removed from patch
120 scalar dMassPatch = 0;
121
122 for (const label celli : pp.faceCells())
123 {
124 scalar ddelta = max(0.0, delta[celli] - deltaStable_);
125 scalar dMass = ddelta*rho[celli]*magSf[celli];
126 massToInject[celli] += dMass;
127 availableMass[celli] -= dMass;
128 dMassPatch += dMass;
129 }
130
131 patchInjectedMasses_[pidi] += dMassPatch;
132 addToInjectedMass(dMassPatch);
133 }
134
136
137 if (writeTime())
138 {
139 scalarField patchInjectedMasses0
140 (
142 (
143 "patchInjectedMasses",
145 )
146 );
147
150 patchInjectedMasses0 += patchInjectedMassTotals;
151
153 (
154 "patchInjectedMasses",
155 patchInjectedMasses0
156 );
157
159 }
160}
161
162
164{
165 // Do not correct if no patches selected
166 if (!patchIDs_.size()) return;
167
168 scalarField patchInjectedMasses
169 (
171 (
172 "patchInjectedMasses",
174 )
175 );
176
179
180 forAll(patchIDs_, pidi)
181 {
182 label patchi = patchIDs_[pidi];
183 patchMasses[patchi] +=
184 patchInjectedMasses[pidi] + patchInjectedMassTotals[pidi];
185 }
186}
187
188
189// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190
191} // End namespace surfaceFilmModels
192} // End namespace regionModels
193} // End namespace Foam
194
195// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
static void listCombineGather(UList< T > &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::listGather with an in-place cop.
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const fvMesh & regionMesh() const
Return the region mesh database.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
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.
Base class for film injection models, handling mass transfer from the film.
void addToInjectedMass(const scalar dMass)
Add to injected mass.
Remove and inject the mass in the film as it passes over the selected patches.
patchInjection(const patchInjection &)=delete
No copy construct.
scalarField patchInjectedMasses_
Injected mass for each patch at which the film is removed.
scalar deltaStable_
Stable film thickness - mass only removed if thickness exceeds.
virtual void patchInjectedMassTotals(scalarField &patchMasses) const
Accumulate the total mass injected for the patches into the.
labelList patchIDs_
List of patch IDs at which the film is removed.
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
label nPatches
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299