pressureDirectedInletVelocityFvPatchVectorField.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-2016 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
38Foam::pressureDirectedInletVelocityFvPatchVectorField::
39pressureDirectedInletVelocityFvPatchVectorField
40(
41 const fvPatch& p,
43)
44:
45 fixedValueFvPatchVectorField(p, iF),
46 phiName_("phi"),
47 rhoName_("rho"),
48 inletDir_(p.size())
49{}
50
51
52Foam::pressureDirectedInletVelocityFvPatchVectorField::
53pressureDirectedInletVelocityFvPatchVectorField
54(
56 const fvPatch& p,
58 const fvPatchFieldMapper& mapper
59)
60:
61 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
62 phiName_(ptf.phiName_),
63 rhoName_(ptf.rhoName_),
64 inletDir_(ptf.inletDir_, mapper)
65{}
66
67
68Foam::pressureDirectedInletVelocityFvPatchVectorField::
69pressureDirectedInletVelocityFvPatchVectorField
70(
71 const fvPatch& p,
73 const dictionary& dict
74)
75:
76 fixedValueFvPatchVectorField(p, iF, dict),
77 phiName_(dict.getOrDefault<word>("phi", "phi")),
78 rhoName_(dict.getOrDefault<word>("rho", "rho")),
79 inletDir_("inletDirection", dict, p.size())
80{}
81
82
83Foam::pressureDirectedInletVelocityFvPatchVectorField::
84pressureDirectedInletVelocityFvPatchVectorField
85(
87)
88:
89 fixedValueFvPatchVectorField(pivpvf),
90 phiName_(pivpvf.phiName_),
91 rhoName_(pivpvf.rhoName_),
92 inletDir_(pivpvf.inletDir_)
93{}
94
95
96Foam::pressureDirectedInletVelocityFvPatchVectorField::
97pressureDirectedInletVelocityFvPatchVectorField
98(
101)
102:
103 fixedValueFvPatchVectorField(pivpvf, iF),
104 phiName_(pivpvf.phiName_),
105 rhoName_(pivpvf.rhoName_),
106 inletDir_(pivpvf.inletDir_)
107{}
108
109
110// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111
113(
114 const fvPatchFieldMapper& m
115)
116{
117 fixedValueFvPatchVectorField::autoMap(m);
118 inletDir_.autoMap(m);
119}
120
121
123(
124 const fvPatchVectorField& ptf,
125 const labelList& addr
126)
127{
128 fixedValueFvPatchVectorField::rmap(ptf, addr);
129
131 refCast<const pressureDirectedInletVelocityFvPatchVectorField>(ptf);
132
133 inletDir_.rmap(tiptf.inletDir_, addr);
134}
135
136
138{
139 if (updated())
140 {
141 return;
142 }
143
144 const surfaceScalarField& phi =
145 db().lookupObject<surfaceScalarField>(phiName_);
146
147 const fvsPatchField<scalar>& phip =
148 patch().patchField<surfaceScalarField, scalar>(phi);
149
150 tmp<vectorField> n = patch().nf();
151 tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
152
154 {
155 operator==(inletDir_*phip/ndmagS);
156 }
158 {
159 const fvPatchField<scalar>& rhop =
160 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
161
162 operator==(inletDir_*phip/(rhop*ndmagS));
163 }
164 else
165 {
167 << "dimensions of phi are not correct"
168 << "\n on patch " << this->patch().name()
169 << " of field " << this->internalField().name()
170 << " in file " << this->internalField().objectPath()
171 << exit(FatalError);
172 }
173
174 fixedValueFvPatchVectorField::updateCoeffs();
175}
176
177
179(
180 Ostream& os
181) const
182{
184 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
185 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
186 inletDir_.writeEntry("inletDirection", os);
187 writeEntry("value", os);
188}
189
190
191// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
192
194(
195 const fvPatchField<vector>& pvf
196)
197{
198 fvPatchField<vector>::operator=(inletDir_*(inletDir_ & pvf));
199}
200
201
202// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203
204namespace Foam
205{
207 (
210 );
211}
212
213// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
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:251
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
This velocity inlet boundary condition is applied to patches where the pressure is specified....
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
Foam::surfaceFields.