pressureDirectedInletVelocityFvPatchVectorField.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 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchVectorField(p, iF),
43  phiName_("phi"),
44  rhoName_("rho"),
45  inletDir_(p.size())
46 {}
47 
48 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
57  fixedValueFvPatchVectorField(p, iF, dict),
58  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
59  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
60  inletDir_("inletDirection", dict, p.size())
61 {}
62 
63 
66 (
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
73  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
74  phiName_(ptf.phiName_),
75  rhoName_(ptf.rhoName_),
76  inletDir_(mapper(ptf.inletDir_))
77 {}
78 
79 
82 (
85 )
86 :
87  fixedValueFvPatchVectorField(pivpvf, iF),
88  phiName_(pivpvf.phiName_),
89  rhoName_(pivpvf.rhoName_),
90  inletDir_(pivpvf.inletDir_)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 (
98  const fvPatchFieldMapper& m
99 )
100 {
101  fixedValueFvPatchVectorField::autoMap(m);
102  m(inletDir_, inletDir_);
103 }
104 
105 
107 (
108  const fvPatchVectorField& ptf,
109  const labelList& addr
110 )
111 {
112  fixedValueFvPatchVectorField::rmap(ptf, addr);
113 
115  refCast<const pressureDirectedInletVelocityFvPatchVectorField>(ptf);
116 
117  inletDir_.rmap(tiptf.inletDir_, addr);
118 }
119 
120 
122 (
123  const fvPatchVectorField& ptf
124 )
125 {
126  fixedValueFvPatchVectorField::reset(ptf);
127 
129  refCast<const pressureDirectedInletVelocityFvPatchVectorField>(ptf);
130 
131  inletDir_.reset(tiptf.inletDir_);
132 }
133 
134 
136 {
137  if (updated())
138  {
139  return;
140  }
141 
142  const surfaceScalarField& phi =
143  db().lookupObject<surfaceScalarField>(phiName_);
144 
145  const fvsPatchField<scalar>& phip =
146  patch().patchField<surfaceScalarField, scalar>(phi);
147 
148  tmp<vectorField> n = patch().nf();
149  tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
150 
151  if (phi.dimensions() == dimFlux)
152  {
153  operator==(inletDir_*phip/ndmagS);
154  }
155  else if (phi.dimensions() == dimMassFlux)
156  {
157  const fvPatchField<scalar>& rhop =
158  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
159 
160  operator==(inletDir_*phip/(rhop*ndmagS));
161  }
162  else
163  {
165  << "dimensions of phi are not correct"
166  << "\n on patch " << this->patch().name()
167  << " of field " << this->internalField().name()
168  << " in file " << this->internalField().objectPath()
169  << exit(FatalError);
170  }
171 
172  fixedValueFvPatchVectorField::updateCoeffs();
173 }
174 
175 
177 (
178  Ostream& os
179 ) const
180 {
182  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
183  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
184  writeEntry(os, "inletDirection", inletDir_);
185  writeEntry(os, "value", *this);
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
190 
191 void Foam::pressureDirectedInletVelocityFvPatchVectorField::operator=
192 (
193  const fvPatchField<vector>& pvf
194 )
195 {
196  fvPatchField<vector>::operator=(inletDir_*(inletDir_ & pvf));
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 namespace Foam
203 {
205  (
208  );
209 }
210 
211 // ************************************************************************* //
Foam::surfaceFields.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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
pressureDirectedInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet & dimensions() const
Return dimensions.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
phi
Definition: correctPhi.H:3
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
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimMassFlux
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:258
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This velocity inlet boundary condition is applied to patches where the pressure is specified...
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.