Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
38(
39 const fvPatch& p,
41)
42:
43 mixedFvPatchVectorField(p, iF),
44 phiName_("phi"),
45 rhoName_("rho"),
46 alphaName_("none"),
47 alphaMin_(1.0)
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{
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
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 }
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
195void 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 );
208
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212namespace Foam
213{
215 (
217 pressurePermeableAlphaInletOutletVelocityFvPatchVectorField
218 );
219}
220
221// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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:346
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void readDict(const dictionary &dict)
Read dictionary entries.
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...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual void write(Ostream &) const
Write.
The pressurePermeableAlphaInletOutletVelocity is a velocity inlet-outlet boundary condition which can...
pressurePermeableAlphaInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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,...
const std::string patch
OpenFOAM patch number as a std::string.
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:40
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
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...
const dimensionSet dimVolume(pow3(dimLength))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< vector > fvPatchVectorField
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.