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-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
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 
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 
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 {
87  this->readValueEntry(dict, IOobjectOption::MUST_READ);
88  refValue() = Zero;
89  refGrad() = Zero;
90  valueFraction() = 1.0;
91 }
92 
93 
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 
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 
127 updateCoeffs()
128 {
129  if (updated())
130  {
131  return;
132  }
133 
134  const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
135 
136  const vectorField n(patch().nf());
137 
138  if (phip.internalField().dimensions() == dimVolume/dimTime)
139  {
140  refValue() = (phip/patch().magSf())*n;
141  }
142  else if (phip.internalField().dimensions() == dimMass/dimTime)
143  {
144  const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName_);
145 
146  refValue() = (phip/(rhop*patch().magSf()))*n;
147  }
148  else
149  {
151  << "dimensions of phi are not correct"
152  << "\n on patch " << this->patch().name()
153  << " of field " << this->internalField().name()
154  << " in file " << this->internalField().objectPath()
155  << exit(FatalError);
156  }
157 
158  valueFraction() = neg(phip);
159 
160  if (alphaName_ != "none")
161  {
162  const scalarField& alphap =
163  patch().lookupPatchField<volScalarField>(alphaName_);
164 
165  const scalarField alphaCut(pos(alphap - alphaMin_));
166  valueFraction() = max(alphaCut, valueFraction());
167  forAll (*this, faceI)
168  {
169  if (valueFraction()[faceI] == 1.0)
170  {
171  refValue()[faceI] = Zero;
172  }
173  }
174  }
175 
176  mixedFvPatchVectorField::updateCoeffs();
177 }
178 
179 
181 (
182  Ostream& os
183 ) const
184 {
186  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
187  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
188  os.writeEntryIfDifferent<word>("alpha", "none", alphaName_);
189  os.writeEntryIfDifferent<scalar>("alphaMin", 1, alphaMin_);
190 }
191 
192 
193 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
194 
195 void Foam::pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
196 ::operator=
197 (
198  const fvPatchField<vector>& pvf
199 )
200 {
201  tmp<vectorField> n = patch().nf();
202 
204  (
205  lerp(pvf, n()*(n() & pvf), valueFraction())
206  );
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 namespace Foam
213 {
215  (
217  pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
218  );
219 }
220 
221 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void readDict(const dictionary &dict)
Read dictionary entries.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
The pressurePermeableAlphaInletOutletVelocity is a velocity inlet-outlet boundary condition which can...
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
pressurePermeableAlphaInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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:336
OBJstream os(runTime.globalPath()/outputName)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void write(Ostream &) const
Write.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127