Loading...
Searching...
No Matches
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
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)
51 refValue() = 0.0;
52 refGrad() = 0.0;
53 valueFraction() = 0.0;
54}
55
56
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 (!this->readValueEntry(dict))
79 {
81 }
82}
83
84
87(
88 const prghPermeableAlphaTotalPressureFvPatchScalarField& ptf,
89 const fvPatch& p,
90 const DimensionedField<scalar, volMesh>& iF,
91 const fvPatchFieldMapper& mapper
92)
93:
94 mixedFvPatchField<scalar>(ptf, p, iF, mapper),
95 p0_(ptf.p0_.clone(p.patch())),
96 phiName_(ptf.phiName_),
97 rhoName_(ptf.rhoName_),
98 UName_(ptf.UName_),
99 alphaName_(ptf.alphaName_),
100 alphaMin_(ptf.alphaMin_),
101 curTimeIndex_(-1)
102{}
103
104
107(
109)
110:
111 mixedFvPatchField<scalar>(tppsf),
112 p0_(tppsf.p0_.clone(this->patch().patch())),
113 phiName_(tppsf.phiName_),
114 rhoName_(tppsf.rhoName_),
115 UName_(tppsf.UName_),
116 alphaName_(tppsf.alphaName_),
117 alphaMin_(tppsf.alphaMin_),
118 curTimeIndex_(-1)
119{}
120
121
124(
127)
128:
129 mixedFvPatchField<scalar>(tppsf, iF),
130 p0_(tppsf.p0_.clone(this->patch().patch())),
131 phiName_(tppsf.phiName_),
132 rhoName_(tppsf.rhoName_),
133 UName_(tppsf.UName_),
134 alphaName_(tppsf.alphaName_),
135 alphaMin_(tppsf.alphaMin_),
136 curTimeIndex_(-1)
137{}
138
139
140// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141
143(
144 const fvPatchFieldMapper& m
145)
146{
148
149 if (p0_)
150 {
151 p0_->autoMap(m);
152 }
153}
154
155
157(
158 const fvPatchScalarField& ptf,
159 const labelList& addr
160)
161{
163
164 const auto& tptf =
166
167 if (p0_)
168 {
169 p0_->rmap(tptf.p0_(), addr);
170 }
171}
172
173
175(
176 const scalarField& snGradp
177)
178{
179 if (updated())
180 {
181 return;
182 }
183
184 const scalarField& rhop =
185 patch().lookupPatchField<volScalarField>(rhoName_);
186
187 const scalarField& phip =
188 patch().lookupPatchField<surfaceScalarField>(phiName_);
189
190 const vectorField& Up =
191 patch().lookupPatchField<volVectorField>(UName_);
192
194 meshObjects::gravity::New(db().time());
195
196 const auto& hRef =
197 db().lookupObject<uniformDimensionedScalarField>("hRef");
198
199 const dimensionedScalar ghRef
200 (
201 mag(g.value()) > SMALL
202 ? g & (cmptMag(g.value())/mag(g.value()))*hRef
204 );
205
206 const scalar t = db().time().timeOutputValue();
207
209 (
210 p0_->value(t)
211 - 0.5*rhop*(neg(phip))*magSqr(Up)
212 - rhop*((g.value() & patch().Cf()) - ghRef.value())
213 );
214
215 refValue() = p;
216
217 refGrad() = snGradp;
218
219 if (alphaName_ != "none")
220 {
221 const scalarField& alphap =
222 patch().lookupPatchField<volScalarField>(alphaName_);
223
224 tmp<scalarField> alphaCut(pos(alphap - alphaMin_));
225 valueFraction() = 1 - alphaCut;
226 }
227
228 if (debug)
229 {
230 const scalar phi = gSum(-phip);
231 Info<< valueFraction() << endl;
232 Info<< patch().boundaryMesh().mesh().name() << ':'
233 << patch().name() << ':'
234 << this->internalField().name() << " :"
235 << " mass flux[Kg/s]:" << phi
236 << endl;
237 }
239 curTimeIndex_ = this->db().time().timeIndex();
240
242}
243
244
246{
247 if (updated())
248 {
249 return;
250 }
251
252 if (curTimeIndex_ != this->db().time().timeIndex())
253 {
255 << "updateCoeffs(const scalarField& snGradp) MUST be called before"
256 " updateCoeffs() or evaluate() to set the boundary gradient."
257 << exit(FatalError);
258 }
259}
260
261
263(
264 Ostream& os
265) const
266{
268 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
269 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
270 os.writeEntryIfDifferent<word>("U", "U", UName_);
271 os.writeEntryIfDifferent<word>("alpha", "none", alphaName_);
272 os.writeEntryIfDifferent<scalar>("alphaMin", 1, alphaMin_);
273
274 if (p0_)
275 {
276 p0_->writeData(os);
277 }
279
280
281// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282
283namespace Foam
284{
286 (
288 prghPermeableAlphaTotalPressureFvPatchScalarField
289 );
290
291}
292
293// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const uniformDimensionedVectorField & g
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
const Type & value() const noexcept
Return const reference to value.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
static tmp< fvPatchField< scalar > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< scalar, volMesh > &)
const DimensionedField< scalar, volMesh > & internalField() const noexcept
virtual void operator=(const UList< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
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.
mixedFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
virtual void write(Ostream &) const
Write.
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.
prghPermeableAlphaTotalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal 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:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< scalar > fvPatchScalarField
label timeIndex
dictionary dict