Loading...
Searching...
No Matches
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
32#include "volFields.H"
33#include "surfaceFields.H"
34
35
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
40(
41 const fvPatch& p,
43)
44:
45 fixedValueFvPatchVectorField(p, iF),
46 phiName_("phi"),
47 rhoName_("rho"),
48 inletDir_(p.size())
49{}
50
51
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
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
85(
87)
88:
89 fixedValueFvPatchVectorField(pivpvf),
90 phiName_(pivpvf.phiName_),
91 rhoName_(pivpvf.rhoName_),
92 inletDir_(pivpvf.inletDir_)
93{}
94
95
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
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
132
133 inletDir_.rmap(tiptf.inletDir_, addr);
134}
135
136
138{
139 if (updated())
140 {
141 return;
142 }
143
144 const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
145
146 tmp<vectorField> n = patch().nf();
147 tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
148
149 if (phip.internalField().dimensions() == dimVolume/dimTime)
150 {
151 operator==(inletDir_*phip/ndmagS);
152 }
153 else if (phip.internalField().dimensions() == dimMass/dimTime)
154 {
155 const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName_);
156
157 operator==(inletDir_*phip/(rhop*ndmagS));
158 }
159 else
160 {
162 << "dimensions of phi are not correct"
163 << "\n on patch " << this->patch().name()
164 << " of field " << this->internalField().name()
165 << " in file " << this->internalField().objectPath()
167 }
168
169 fixedValueFvPatchVectorField::updateCoeffs();
170}
171
172
174(
175 Ostream& os
176) const
177{
179 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
180 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
181 inletDir_.writeEntry("inletDirection", os);
183}
184
185
186// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
187
188void Foam::pressureDirectedInletVelocityFvPatchVectorField::operator=
189(
190 const fvPatchField<vector>& pvf
191)
192{
193 fvPatchField<vector>::operator=(inletDir_*(inletDir_ & pvf));
195
196
197// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198
199namespace Foam
200{
202 (
204 pressureDirectedInletVelocityFvPatchVectorField
205 );
206}
207
208// ************************************************************************* //
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
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.
pressureDirectedInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< vector > fvPatchVectorField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict
Foam::surfaceFields.