prghPermeableAlphaTotalPressureFvPatchScalarField.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) 2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
30#include "fvPatchFieldMapper.H"
31#include "gravityMeshObject.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35Foam::prghPermeableAlphaTotalPressureFvPatchScalarField::
36prghPermeableAlphaTotalPressureFvPatchScalarField
37(
38 const fvPatch& p,
40)
41:
42 mixedFvPatchField<scalar>(p, iF),
43 p0_(nullptr),
44 phiName_("phi"),
45 rhoName_("rho"),
46 UName_("U"),
47 alphaName_("none"),
48 alphaMin_(1.0),
49 curTimeIndex_(-1)
50{
51 refValue() = 0.0;
52 refGrad() = 0.0;
53 valueFraction() = 0.0;
54}
55
56
57Foam::prghPermeableAlphaTotalPressureFvPatchScalarField::
58prghPermeableAlphaTotalPressureFvPatchScalarField
59(
60 const fvPatch& p,
62 const dictionary& dict
63)
64:
65 mixedFvPatchField<scalar>(p, iF),
66 p0_(PatchFunction1<scalar>::New(p.patch(), "p", dict)),
67 phiName_(dict.getOrDefault<word>("phi", "phi")),
68 rhoName_(dict.getOrDefault<word>("rho", "rho")),
69 UName_(dict.getOrDefault<word>("U", "U")),
70 alphaName_(dict.getOrDefault<word>("alpha", "none")),
71 alphaMin_(dict.getOrDefault<scalar>("alphaMin", 1)),
72 curTimeIndex_(-1)
73{
74 refValue() = 1.0;
75 refGrad() = 0.0;
76 valueFraction() = 0.0;
77
78 if (dict.found("value"))
79 {
81 (
82 Field<scalar>("value", dict, p.size())
83 );
84 }
85 else
86 {
88 }
89}
90
91
92Foam::prghPermeableAlphaTotalPressureFvPatchScalarField::
93prghPermeableAlphaTotalPressureFvPatchScalarField
94(
96 const fvPatch& p,
98 const fvPatchFieldMapper& mapper
99)
100:
101 mixedFvPatchField<scalar>(ptf, p, iF, mapper),
102 p0_(ptf.p0_.clone(p.patch())),
103 phiName_(ptf.phiName_),
104 rhoName_(ptf.rhoName_),
105 UName_(ptf.UName_),
106 alphaName_(ptf.alphaName_),
107 alphaMin_(ptf.alphaMin_),
108 curTimeIndex_(-1)
109{}
110
111
112Foam::prghPermeableAlphaTotalPressureFvPatchScalarField::
113prghPermeableAlphaTotalPressureFvPatchScalarField
114(
116)
117:
118 mixedFvPatchField<scalar>(tppsf),
119 p0_(tppsf.p0_.clone(this->patch().patch())),
120 phiName_(tppsf.phiName_),
121 rhoName_(tppsf.rhoName_),
122 UName_(tppsf.UName_),
123 alphaName_(tppsf.alphaName_),
124 alphaMin_(tppsf.alphaMin_),
125 curTimeIndex_(-1)
126{}
127
128
129Foam::prghPermeableAlphaTotalPressureFvPatchScalarField::
130prghPermeableAlphaTotalPressureFvPatchScalarField
131(
134)
135:
136 mixedFvPatchField<scalar>(tppsf, iF),
137 p0_(tppsf.p0_.clone(this->patch().patch())),
138 phiName_(tppsf.phiName_),
139 rhoName_(tppsf.rhoName_),
140 UName_(tppsf.UName_),
141 alphaName_(tppsf.alphaName_),
142 alphaMin_(tppsf.alphaMin_),
143 curTimeIndex_(-1)
144{}
145
146
147// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148
150(
151 const fvPatchFieldMapper& m
152)
153{
155
156 if (p0_)
157 {
158 p0_->autoMap(m);
159 }
160}
161
162
164(
165 const fvPatchScalarField& ptf,
166 const labelList& addr
167)
168{
170
171 const auto& tptf =
172 refCast<const prghPermeableAlphaTotalPressureFvPatchScalarField>(ptf);
173
174 if (p0_)
175 {
176 p0_->rmap(tptf.p0_(), addr);
177 }
178}
179
180
182(
183 const scalarField& snGradp
184)
185{
186 if (updated())
187 {
188 return;
189 }
190
191 const scalarField& rhop =
192 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
193
194 const scalarField& phip =
195 patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
196
197 const vectorField& Up =
198 patch().lookupPatchField<volVectorField, vector>(UName_);
199
201 meshObjects::gravity::New(db().time());
202
203 const auto& hRef =
204 db().lookupObject<uniformDimensionedScalarField>("hRef");
205
206 const dimensionedScalar ghRef
207 (
208 mag(g.value()) > SMALL
209 ? g & (cmptMag(g.value())/mag(g.value()))*hRef
211 );
212
213 const scalar t = db().time().timeOutputValue();
214
216 (
217 p0_->value(t)
218 - 0.5*rhop*(1.0 - pos0(phip))*magSqr(Up)
219 - rhop*((g.value() & patch().Cf()) - ghRef.value())
220 );
221
222 refValue() = p;
223
224 refGrad() = snGradp;
225
226 if (alphaName_ != "none")
227 {
228 const scalarField& alphap =
229 patch().lookupPatchField<volScalarField, scalar>(alphaName_);
230 tmp<scalarField> alphaCut(pos(alphap - alphaMin_));
231 valueFraction() = 1 - alphaCut;
232 }
233
234 if (debug)
235 {
236 const scalar phi = gSum(-phip);
237 Info<< valueFraction() << endl;
238 Info<< patch().boundaryMesh().mesh().name() << ':'
239 << patch().name() << ':'
240 << this->internalField().name() << " :"
241 << " mass flux[Kg/s]:" << phi
242 << endl;
243 }
244
245 curTimeIndex_ = this->db().time().timeIndex();
246
248}
249
250
252{
253 if (updated())
254 {
255 return;
256 }
257
258 if (curTimeIndex_ != this->db().time().timeIndex())
259 {
261 << "updateCoeffs(const scalarField& snGradp) MUST be called before"
262 " updateCoeffs() or evaluate() to set the boundary gradient."
263 << exit(FatalError);
264 }
265}
266
267
269(
270 Ostream& os
271) const
272{
274 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
275 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
276 os.writeEntryIfDifferent<word>("U", "U", UName_);
277 os.writeEntryIfDifferent<word>("alpha", "none", alphaName_);
278 os.writeEntryIfDifferent<scalar>("alphaMin", 1, alphaMin_);
279
280 if (p0_)
281 {
282 p0_->writeData(os);
283 }
284}
285
286
287// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288
289namespace Foam
290{
292 (
295 );
296
297}
298
299// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const uniformDimensionedVectorField & g
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type.
Definition: Field.H:82
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual Field< scalar > & refGrad()
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual Field< scalar > & refValue()
virtual scalarField & valueFraction()
The prghPermeableAlphaTotalPressure is a mixed boundary condition for the p_rgh variable in multiphas...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void updateSnGrad(const scalarField &snGradp)
Update the patch pressure gradient field from the given snGradp.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict