fixedNormalInletOutletVelocityFvPatchVectorField.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) 2014-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
37Foam::fixedNormalInletOutletVelocityFvPatchVectorField::
38fixedNormalInletOutletVelocityFvPatchVectorField
39(
40 const fvPatch& p,
42)
43:
44 directionMixedFvPatchVectorField(p, iF),
45 phiName_("phi"),
46 fixTangentialInflow_(true),
47 normalVelocity_
48 (
49 fvPatchVectorField::New("fixedValue", p, iF)
50 )
51{
52 refValue() = Zero;
53 refGrad() = Zero;
54 valueFraction() = Zero;
55}
56
57
58Foam::fixedNormalInletOutletVelocityFvPatchVectorField::
59fixedNormalInletOutletVelocityFvPatchVectorField
60(
61 const fvPatch& p,
63 const dictionary& dict
64)
65:
66 directionMixedFvPatchVectorField(p, iF),
67 phiName_(dict.getOrDefault<word>("phi", "phi")),
68 fixTangentialInflow_(dict.lookup("fixTangentialInflow")),
69 normalVelocity_
70 (
71 fvPatchVectorField::New(p, iF, dict.subDict("normalVelocity"))
72 )
73{
74 patchType() = dict.getOrDefault<word>("patchType", word::null);
76 refValue() = normalVelocity();
77 refGrad() = Zero;
78 valueFraction() = Zero;
79}
80
81
82Foam::fixedNormalInletOutletVelocityFvPatchVectorField::
83fixedNormalInletOutletVelocityFvPatchVectorField
84(
86 const fvPatch& p,
88 const fvPatchFieldMapper& mapper
89)
90:
91 directionMixedFvPatchVectorField(ptf, p, iF, mapper),
92 phiName_(ptf.phiName_),
93 fixTangentialInflow_(ptf.fixTangentialInflow_),
94 normalVelocity_
95 (
96 fvPatchVectorField::New(ptf.normalVelocity(), p, iF, mapper)
97 )
98{}
99
100
101Foam::fixedNormalInletOutletVelocityFvPatchVectorField::
102fixedNormalInletOutletVelocityFvPatchVectorField
103(
105)
106:
107 directionMixedFvPatchVectorField(pivpvf),
108 phiName_(pivpvf.phiName_),
109 fixTangentialInflow_(pivpvf.fixTangentialInflow_),
110 normalVelocity_(pivpvf.normalVelocity().clone())
111{}
112
113
114Foam::fixedNormalInletOutletVelocityFvPatchVectorField::
115fixedNormalInletOutletVelocityFvPatchVectorField
116(
119)
120:
121 directionMixedFvPatchVectorField(pivpvf, iF),
122 phiName_(pivpvf.phiName_),
123 fixTangentialInflow_(pivpvf.fixTangentialInflow_),
124 normalVelocity_(pivpvf.normalVelocity().clone())
125{}
126
127
128// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129
131(
132 const fvPatchFieldMapper& m
133)
134{
135 directionMixedFvPatchVectorField::autoMap(m);
136 normalVelocity_->autoMap(m);
137}
138
139
141(
142 const fvPatchVectorField& ptf,
143 const labelList& addr
144)
145{
146 directionMixedFvPatchVectorField::rmap(ptf, addr);
147
149 refCast<const fixedNormalInletOutletVelocityFvPatchVectorField>(ptf);
150
151 normalVelocity_->rmap(fniovptf.normalVelocity(), addr);
152}
153
154
156{
157 if (updated())
158 {
159 return;
160 }
161
162 normalVelocity_->evaluate();
163 refValue() = normalVelocity();
164
165 valueFraction() = sqr(patch().nf());
166
167 if (fixTangentialInflow_)
168 {
169 const fvsPatchField<scalar>& phip =
170 patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
171
172 valueFraction() += neg(phip)*(I - valueFraction());
173 }
174
175 directionMixedFvPatchVectorField::updateCoeffs();
176 directionMixedFvPatchVectorField::evaluate();
177}
178
179
181(
182 Ostream& os
183)
184const
185{
187 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
188 os.writeEntry("fixTangentialInflow", fixTangentialInflow_);
189
190 os.beginBlock("normalVelocity");
191 normalVelocity_->write(os);
192 os.endBlock();
193
194 writeEntry("value", os);
195}
196
197
198// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
199
201(
202 const fvPatchField<vector>& pvf
203)
204{
205 tmp<vectorField> normalValue = transform(valueFraction(), refValue());
206 tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
207 fvPatchField<vector>::operator=(normalValue + transformGradValue);
208}
209
210
211// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212
213namespace Foam
214{
216 (
219 );
220}
221
222// ************************************************************************* //
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:62
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:105
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:87
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
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
This velocity inlet/outlet boundary condition combines a fixed normal component obtained from the "no...
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.
const fvPatchVectorField & normalVelocity() const
Return the BC which provides the normal component of velocity.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< vector > &)
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
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const Identity< scalar > I
Definition: Identity.H:94
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:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
Foam::surfaceFields.