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-2024 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 "fieldMapper.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", dimless, dict, p.size())
46 {
48  (
49  vectorField("value", iF.dimensions(), dict, p.size())
50  );
51  refValue() = *this;
52  refGrad() = Zero;
53  valueFraction() = 0.0;
54 }
55 
56 
59 (
61  const fvPatch& p,
63  const fieldMapper& mapper
64 )
65 :
66  mixedFvPatchVectorField(ptf, p, iF, mapper),
67  phiName_(ptf.phiName_),
68  rhoName_(ptf.rhoName_),
69  inletDir_(mapper(ptf.inletDir_))
70 {}
71 
72 
75 (
78 )
79 :
80  mixedFvPatchVectorField(pivpvf, iF),
81  phiName_(pivpvf.phiName_),
82  rhoName_(pivpvf.rhoName_),
83  inletDir_(pivpvf.inletDir_)
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
90 (
91  const fvPatchVectorField& ptf,
92  const fieldMapper& mapper
93 )
94 {
95  mixedFvPatchVectorField::map(ptf, mapper);
96 
98  refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
99  (ptf);
100 
101  mapper(inletDir_, tiptf.inletDir_);
102 }
103 
104 
106 (
107  const fvPatchVectorField& ptf
108 )
109 {
110  mixedFvPatchVectorField::reset(ptf);
111 
113  refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
114  (ptf);
115 
116  inletDir_.reset(tiptf.inletDir_);
117 }
118 
119 
121 {
122  if (updated())
123  {
124  return;
125  }
126 
127  const surfaceScalarField& phi =
128  db().lookupObject<surfaceScalarField>(phiName_);
129 
130  const fvsPatchField<scalar>& phip =
131  patch().patchField<surfaceScalarField, scalar>(phi);
132 
133  tmp<vectorField> n = patch().nf();
134  tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
135 
136  if (phi.dimensions() == dimVolumetricFlux)
137  {
138  refValue() = inletDir_*phip/ndmagS;
139  }
140  else if (phi.dimensions() == dimMassFlux)
141  {
142  const fvPatchField<scalar>& rhop =
143  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
144 
145  refValue() = inletDir_*phip/(rhop*ndmagS);
146  }
147  else
148  {
150  << "dimensions of phi are not correct"
151  << "\n on patch " << this->patch().name()
152  << " of field " << this->internalField().name()
153  << " in file " << this->internalField().objectPath()
154  << exit(FatalError);
155  }
156 
157  valueFraction() = neg(phip);
158 
159  mixedFvPatchVectorField::updateCoeffs();
160 }
161 
162 
164 (
165  Ostream& os
166 ) const
167 {
169  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
170  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
171  writeEntry(os, "inletDirection", inletDir_);
172  writeEntry(os, "value", *this);
173 }
174 
175 
176 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
177 
178 void Foam::pressureDirectedInletOutletVelocityFvPatchVectorField::operator=
179 (
180  const fvPatchField<vector>& pvf
181 )
182 {
184  (
185  valueFraction()*(inletDir_*(inletDir_ & pvf))
186  + (1 - valueFraction())*pvf
187  );
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 namespace Foam
194 {
196  (
199  );
200 }
201 
202 // ************************************************************************* //
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:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
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:83
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.
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.
virtual void map(const fvPatchVectorField &, const fieldMapper &)
Map the given fvPatchField onto this 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:334
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
const dimensionSet dimVolumetricFlux
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.