pressureDirectedInletOutletVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  mixedFvPatchVectorField(p, iF),
42  phiName_("phi"),
43  rhoName_("rho"),
44  inletDir_(p.size())
45 {
46  refValue() = *this;
47  refGrad() = Zero;
48  valueFraction() = 0.0;
49 }
50 
51 
54 (
55  const fvPatch& p,
57  const dictionary& dict
58 )
59 :
60  mixedFvPatchVectorField(p, iF),
61  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
62  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
63  inletDir_("inletDirection", dict, p.size())
64 {
65  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
66  refValue() = *this;
67  refGrad() = Zero;
68  valueFraction() = 0.0;
69 }
70 
71 
74 (
76  const fvPatch& p,
78  const fvPatchFieldMapper& mapper
79 )
80 :
81  mixedFvPatchVectorField(ptf, p, iF, mapper),
82  phiName_(ptf.phiName_),
83  rhoName_(ptf.rhoName_),
84  inletDir_(mapper(ptf.inletDir_))
85 {}
86 
87 
90 (
93 )
94 :
95  mixedFvPatchVectorField(pivpvf, iF),
96  phiName_(pivpvf.phiName_),
97  rhoName_(pivpvf.rhoName_),
98  inletDir_(pivpvf.inletDir_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 (
106  const fvPatchFieldMapper& m
107 )
108 {
109  mixedFvPatchVectorField::autoMap(m);
110  m(inletDir_, inletDir_);
111 }
112 
113 
115 (
116  const fvPatchVectorField& ptf,
117  const labelList& addr
118 )
119 {
120  mixedFvPatchVectorField::rmap(ptf, addr);
121 
123  refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
124  (ptf);
125 
126  inletDir_.rmap(tiptf.inletDir_, addr);
127 }
128 
129 
131 (
132  const fvPatchVectorField& ptf
133 )
134 {
135  mixedFvPatchVectorField::reset(ptf);
136 
138  refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
139  (ptf);
140 
141  inletDir_.reset(tiptf.inletDir_);
142 }
143 
144 
146 {
147  if (updated())
148  {
149  return;
150  }
151 
152  const surfaceScalarField& phi =
153  db().lookupObject<surfaceScalarField>(phiName_);
154 
155  const fvsPatchField<scalar>& phip =
156  patch().patchField<surfaceScalarField, scalar>(phi);
157 
158  tmp<vectorField> n = patch().nf();
159  tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
160 
161  if (phi.dimensions() == dimFlux)
162  {
163  refValue() = inletDir_*phip/ndmagS;
164  }
165  else if (phi.dimensions() == dimMassFlux)
166  {
167  const fvPatchField<scalar>& rhop =
168  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
169 
170  refValue() = inletDir_*phip/(rhop*ndmagS);
171  }
172  else
173  {
175  << "dimensions of phi are not correct"
176  << "\n on patch " << this->patch().name()
177  << " of field " << this->internalField().name()
178  << " in file " << this->internalField().objectPath()
179  << exit(FatalError);
180  }
181 
182  valueFraction() = neg(phip);
183 
184  mixedFvPatchVectorField::updateCoeffs();
185 }
186 
187 
189 (
190  Ostream& os
191 ) const
192 {
194  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
195  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
196  writeEntry(os, "inletDirection", inletDir_);
197  writeEntry(os, "value", *this);
198 }
199 
200 
201 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
202 
203 void Foam::pressureDirectedInletOutletVelocityFvPatchVectorField::operator=
204 (
205  const fvPatchField<vector>& pvf
206 )
207 {
209  (
210  valueFraction()*(inletDir_*(inletDir_ & pvf))
211  + (1 - valueFraction())*pvf
212  );
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 namespace Foam
219 {
221  (
224  );
225 }
226 
227 // ************************************************************************* //
Foam::surfaceFields.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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:306
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
volVectorField vectorField(fieldObject, mesh)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
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:58
const dimensionSet & dimensions() const
Return dimensions.
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
phi
Definition: correctPhi.H:3
static const zero Zero
Definition: zero.H:97
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimMassFlux
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
This velocity inlet/outlet boundary condition is applied to pressure boundaries where the pressure is...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
label n
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.