pressurePermeableAlphaInletOutletVelocityFvPatchVectorField.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 "volFields.H"
32#include "surfaceFields.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36Foam::pressurePermeableAlphaInletOutletVelocityFvPatchVectorField::
37pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
38(
39 const fvPatch& p,
41)
42:
43 mixedFvPatchVectorField(p, iF),
44 phiName_("phi"),
45 rhoName_("rho"),
46 alphaName_("none"),
47 alphaMin_(1.0)
48{
49 refValue() = Zero;
50 refGrad() = Zero;
51 valueFraction() = 1.0;
52}
53
54
55Foam::pressurePermeableAlphaInletOutletVelocityFvPatchVectorField::
56pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
57(
59 const fvPatch& p,
61 const fvPatchFieldMapper& mapper
62)
63:
64 mixedFvPatchVectorField(ptf, p, iF, mapper),
65 phiName_(ptf.phiName_),
66 rhoName_(ptf.rhoName_),
67 alphaName_(ptf.alphaName_),
68 alphaMin_(ptf.alphaMin_)
69{}
70
71
72Foam::pressurePermeableAlphaInletOutletVelocityFvPatchVectorField::
73pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
74(
75 const fvPatch& p,
77 const dictionary& dict
78)
79:
80 mixedFvPatchVectorField(p, iF),
81 phiName_(dict.getOrDefault<word>("phi", "phi")),
82 rhoName_(dict.getOrDefault<word>("rho", "rho")),
83 alphaName_(dict.getOrDefault<word>("alpha", "none")),
84 alphaMin_(dict.getOrDefault<scalar>("alphaMin", 1))
85{
86 patchType() = dict.getOrDefault<word>("patchType", word::null);
88 refValue() = Zero;
89 refGrad() = Zero;
90 valueFraction() = 1.0;
91}
92
93
94Foam::pressurePermeableAlphaInletOutletVelocityFvPatchVectorField::
95pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
96(
98)
99:
100 mixedFvPatchVectorField(pivpvf),
101 phiName_(pivpvf.phiName_),
102 rhoName_(pivpvf.rhoName_),
103 alphaName_(pivpvf.alphaName_),
104 alphaMin_(pivpvf.alphaMin_)
105{}
106
107
108Foam::pressurePermeableAlphaInletOutletVelocityFvPatchVectorField::
109pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
110(
113)
114:
115 mixedFvPatchVectorField(pivpvf, iF),
116 phiName_(pivpvf.phiName_),
117 rhoName_(pivpvf.rhoName_),
118 alphaName_(pivpvf.alphaName_),
119 alphaMin_(pivpvf.alphaMin_)
120{}
121
122
123// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124
125
128{
129 if (updated())
130 {
131 return;
132 }
133
134 const auto& phi = db().lookupObject<surfaceScalarField>(phiName_);
135
136 const fvsPatchField<scalar>& phip =
137 patch().patchField<surfaceScalarField, scalar>(phi);
138
139 const vectorField n(patch().nf());
140
142 {
143 refValue() = (phip/patch().magSf())*n;
144 }
146 {
147 const fvPatchField<scalar>& rhop =
148 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
149
150 refValue() = (phip/(rhop*patch().magSf()))*n;
151 }
152 else
153 {
155 << "dimensions of phi are not correct"
156 << "\n on patch " << this->patch().name()
157 << " of field " << this->internalField().name()
158 << " in file " << this->internalField().objectPath()
159 << exit(FatalError);
160 }
161
162 valueFraction() = 1.0 - pos0(phip);
163
164 if (alphaName_ != "none")
165 {
166 const scalarField& alphap =
167 patch().lookupPatchField<volScalarField, scalar>(alphaName_);
168
169 const scalarField alphaCut(pos(alphap - alphaMin_));
170 valueFraction() = max(alphaCut, valueFraction());
171 forAll (*this, faceI)
172 {
173 if (valueFraction()[faceI] == 1.0)
174 {
175 refValue()[faceI] = Zero;
176 }
177 }
178 }
179
180 mixedFvPatchVectorField::updateCoeffs();
181}
182
183
185(
186 Ostream& os
187) const
188{
189 mixedFvPatchVectorField::write(os);
190 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
191 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
192 os.writeEntryIfDifferent<word>("alpha", "none", alphaName_);
193 os.writeEntryIfDifferent<scalar>("alphaMin", 1, alphaMin_);
194}
195
196
197// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
198
199void Foam::pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
200::operator=
201(
202 const fvPatchField<vector>& pvf
203)
204{
205 tmp<vectorField> n = patch().nf();
206
208 (
209 valueFraction()*(n()*(n() & pvf))
210 + (1 - valueFraction())*pvf
211 );
212}
213
214
215// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216
217namespace Foam
218{
220 (
223 );
224}
225
226// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
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
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
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< vector > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
The pressurePermeableAlphaInletOutletVelocity is a velocity inlet-outlet boundary condition which can...
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
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)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.