pressureDirectedInletOutletVelocityFvPatchVectorField.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2017-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  mixedFvPatchVectorField(p, iF),
45  phiName_("phi"),
46  rhoName_("rho"),
47  inletDir_(p.size())
48 {
49  refValue() = *this;
50  refGrad() = Zero;
51  valueFraction() = 0.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  inletDir_(ptf.inletDir_, mapper)
68 {}
69 
70 
73 (
74  const fvPatch& p,
76  const dictionary& dict
77 )
78 :
79  mixedFvPatchVectorField(p, iF),
80  phiName_(dict.getOrDefault<word>("phi", "phi")),
81  rhoName_(dict.getOrDefault<word>("rho", "rho")),
82  inletDir_("inletDirection", dict, p.size())
83 {
85  this->readValueEntry(dict, IOobjectOption::MUST_READ);
86  refValue() = *this;
87  refGrad() = Zero;
88  valueFraction() = 0.0;
89 }
90 
91 
94 (
96 )
97 :
98  mixedFvPatchVectorField(pivpvf),
99  phiName_(pivpvf.phiName_),
100  rhoName_(pivpvf.rhoName_),
101  inletDir_(pivpvf.inletDir_)
102 {}
103 
104 
107 (
110 )
111 :
112  mixedFvPatchVectorField(pivpvf, iF),
113  phiName_(pivpvf.phiName_),
114  rhoName_(pivpvf.rhoName_),
115  inletDir_(pivpvf.inletDir_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 (
123  const fvPatchFieldMapper& m
124 )
125 {
126  mixedFvPatchVectorField::autoMap(m);
127  inletDir_.autoMap(m);
128 }
129 
130 
132 (
133  const fvPatchVectorField& ptf,
134  const labelList& addr
135 )
136 {
137  mixedFvPatchVectorField::rmap(ptf, addr);
138 
140  refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
141  (ptf);
142 
143  inletDir_.rmap(tiptf.inletDir_, addr);
144 }
145 
146 
148 {
149  if (updated())
150  {
151  return;
152  }
153 
154  const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
155 
156  tmp<vectorField> n = patch().nf();
157  tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
158 
159  if (phip.internalField().dimensions() == dimVolume/dimTime)
160  {
161  refValue() = inletDir_*phip/ndmagS;
162  }
163  else if (phip.internalField().dimensions() == dimMass/dimTime)
164  {
165  const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName_);
166 
167  refValue() = inletDir_*phip/(rhop*ndmagS);
168  }
169  else
170  {
172  << "dimensions of phi are not correct"
173  << "\n on patch " << this->patch().name()
174  << " of field " << this->internalField().name()
175  << " in file " << this->internalField().objectPath()
176  << exit(FatalError);
177  }
178 
179  valueFraction() = neg(phip);
180 
181  mixedFvPatchVectorField::updateCoeffs();
182 }
183 
184 
186 (
187  Ostream& os
188 ) const
189 {
191  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
192  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
193  inletDir_.writeEntry("inletDirection", os);
195 }
196 
197 
198 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
199 
200 void Foam::pressureDirectedInletOutletVelocityFvPatchVectorField::operator=
201 (
202  const fvPatchField<vector>& pvf
203 )
204 {
206  (
207  lerp(pvf, inletDir_*(inletDir_ & pvf), valueFraction())
208  );
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 namespace Foam
215 {
217  (
219  pressureDirectedInletOutletVelocityFvPatchVectorField
220  );
221 }
222 
223 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
pressureDirectedInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void readDict(const dictionary &dict)
Read dictionary entries.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
This velocity inlet/outlet boundary condition is applied to velocity boundaries where the pressure is...
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
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127