Loading...
Searching...
No Matches
pressureInletOutletVelocityFvPatchVectorField.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) 2017-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 "volFields.H"
32#include "surfaceFields.H"
33
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
39(
40 const fvPatch& p,
42)
43:
44 directionMixedFvPatchVectorField(p, iF),
45 phiName_("phi")
47 refValue() = Zero;
48 refGrad() = Zero;
49 valueFraction() = Zero;
50}
51
52
55(
57 const fvPatch& p,
59 const fvPatchFieldMapper& mapper
60)
61:
62 directionMixedFvPatchVectorField(ptf, p, iF, mapper),
63 phiName_(ptf.phiName_)
64{
65 if (ptf.tangentialVelocity_.size())
66 {
67 tangentialVelocity_ = mapper(ptf.tangentialVelocity_);
68 }
69}
70
71
74(
75 const fvPatch& p,
77 const dictionary& dict
78)
79:
80 directionMixedFvPatchVectorField(p, iF),
81 phiName_(dict.getOrDefault<word>("phi", "phi"))
82{
84 this->readValueEntry(dict, IOobjectOption::MUST_READ);
85
86 if (dict.found("tangentialVelocity"))
87 {
89 (
90 vectorField("tangentialVelocity", dict, p.size())
91 );
92 }
93 else
94 {
95 refValue() = Zero;
106 const pressureInletOutletVelocityFvPatchVectorField& pivpvf
107)
109 directionMixedFvPatchVectorField(pivpvf),
110 phiName_(pivpvf.phiName_),
111 tangentialVelocity_(pivpvf.tangentialVelocity_)
112{}
113
114
117(
120)
121:
122 directionMixedFvPatchVectorField(pivpvf, iF),
123 phiName_(pivpvf.phiName_),
124 tangentialVelocity_(pivpvf.tangentialVelocity_)
125{}
126
127
128// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129
131setTangentialVelocity(const vectorField& tangentialVelocity)
133 tangentialVelocity_ = tangentialVelocity;
134 const vectorField n(patch().nf());
135 refValue() = tangentialVelocity_ - n*(n & tangentialVelocity_);
136}
137
138
140(
141 const fvPatchFieldMapper& m
142)
143{
144 directionMixedFvPatchVectorField::autoMap(m);
145 if (tangentialVelocity_.size())
146 {
147 tangentialVelocity_.autoMap(m);
148 }
149}
150
151
153(
154 const fvPatchVectorField& ptf,
155 const labelList& addr
156)
157{
158 directionMixedFvPatchVectorField::rmap(ptf, addr);
159
160 if (tangentialVelocity_.size())
161 {
162 const pressureInletOutletVelocityFvPatchVectorField& tiptf =
164
165 tangentialVelocity_.rmap(tiptf.tangentialVelocity_, addr);
166 }
167}
168
169
171{
172 if (updated())
173 {
174 return;
175 }
176
177 const fvsPatchField<scalar>& phip =
178 patch().lookupPatchField<surfaceScalarField>(phiName_);
179
180 valueFraction() = neg(phip)*(I - sqr(patch().nf()));
181
182 directionMixedFvPatchVectorField::updateCoeffs();
183 directionMixedFvPatchVectorField::evaluate();
184}
185
186
188(
189 Ostream& os
190)
191const
192{
194 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
195 if (tangentialVelocity_.size())
196 {
197 tangentialVelocity_.writeEntry("tangentialVelocity", os);
200}
201
202
203// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
204
205void Foam::pressureInletOutletVelocityFvPatchVectorField::operator=
206(
207 const fvPatchField<vector>& pvf
208)
209{
210 tmp<vectorField> normalValue = transform(valueFraction(), refValue());
211 tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
212 fvPatchField<vector>::operator=(normalValue + transformGradValue);
214
215
216// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217
218namespace Foam
219{
221 (
223 pressureInletOutletVelocityFvPatchVectorField
224 );
225}
226
227// ************************************************************************* //
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...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void readDict(const dictionary &dict)
Read dictionary entries.
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
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
This velocity inlet/outlet boundary condition is applied to velocity boundaries where the pressure is...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void setTangentialVelocity(const vectorField &tangentialVelocity)
Reset the tangential velocity.
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.
pressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const vectorField & tangentialVelocity() const
Return the tangential velocity.
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
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
static const Identity< scalar > I
Definition Identity.H:100
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< vector > fvPatchVectorField
dictionary dict
Foam::surfaceFields.