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-2023 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  const dictionary& dict
40 )
41 :
42  mixedFvPatchVectorField(p, iF, dict, false),
43  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
44  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
45  inletDir_("inletDirection", dict, p.size())
46 {
48  refValue() = *this;
49  refGrad() = Zero;
50  valueFraction() = 0.0;
51 }
52 
53 
56 (
58  const fvPatch& p,
60  const fvPatchFieldMapper& mapper
61 )
62 :
63  mixedFvPatchVectorField(ptf, p, iF, mapper),
64  phiName_(ptf.phiName_),
65  rhoName_(ptf.rhoName_),
66  inletDir_(mapper(ptf.inletDir_))
67 {}
68 
69 
72 (
75 )
76 :
77  mixedFvPatchVectorField(pivpvf, iF),
78  phiName_(pivpvf.phiName_),
79  rhoName_(pivpvf.rhoName_),
80  inletDir_(pivpvf.inletDir_)
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
87 (
88  const fvPatchVectorField& ptf,
89  const fvPatchFieldMapper& mapper
90 )
91 {
92  mixedFvPatchVectorField::map(ptf, mapper);
93 
95  refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
96  (ptf);
97 
98  mapper(inletDir_, tiptf.inletDir_);
99 }
100 
101 
103 (
104  const fvPatchVectorField& ptf
105 )
106 {
107  mixedFvPatchVectorField::reset(ptf);
108 
110  refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
111  (ptf);
112 
113  inletDir_.reset(tiptf.inletDir_);
114 }
115 
116 
118 {
119  if (updated())
120  {
121  return;
122  }
123 
124  const surfaceScalarField& phi =
125  db().lookupObject<surfaceScalarField>(phiName_);
126 
127  const fvsPatchField<scalar>& phip =
128  patch().patchField<surfaceScalarField, scalar>(phi);
129 
130  tmp<vectorField> n = patch().nf();
131  tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
132 
133  if (phi.dimensions() == dimFlux)
134  {
135  refValue() = inletDir_*phip/ndmagS;
136  }
137  else if (phi.dimensions() == dimMassFlux)
138  {
139  const fvPatchField<scalar>& rhop =
140  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
141 
142  refValue() = inletDir_*phip/(rhop*ndmagS);
143  }
144  else
145  {
147  << "dimensions of phi are not correct"
148  << "\n on patch " << this->patch().name()
149  << " of field " << this->internalField().name()
150  << " in file " << this->internalField().objectPath()
151  << exit(FatalError);
152  }
153 
154  valueFraction() = neg(phip);
155 
156  mixedFvPatchVectorField::updateCoeffs();
157 }
158 
159 
161 (
162  Ostream& os
163 ) const
164 {
166  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
167  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
168  writeEntry(os, "inletDirection", inletDir_);
169  writeEntry(os, "value", *this);
170 }
171 
172 
173 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
174 
175 void Foam::pressureDirectedInletOutletVelocityFvPatchVectorField::operator=
176 (
177  const fvPatchField<vector>& pvf
178 )
179 {
181  (
182  valueFraction()*(inletDir_*(inletDir_ & pvf))
183  + (1 - valueFraction())*pvf
184  );
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 namespace Foam
191 {
193  (
196  );
197 }
198 
199 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
This velocity inlet/outlet boundary condition is applied to pressure boundaries where the pressure is...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchVectorField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
pressureDirectedInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const dimensionSet dimMassFlux
SurfaceField< scalar > surfaceScalarField
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const dimensionSet dimFlux
dictionary dict
volScalarField & p
Foam::surfaceFields.