Loading...
Searching...
No Matches
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
39(
40 const fvPatch& p,
42)
43:
44 directionMixedFvPatchVectorField(p, iF),
45 phiName_("phi"),
48 (
49 fvPatchVectorField::New("fixedValue", p, iF)
50 )
52 refValue() = Zero;
53 refGrad() = Zero;
54 valueFraction() = Zero;
55}
56
57
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{
75 this->readValueEntry(dict, IOobjectOption::MUST_READ);
76 refValue() = normalVelocity();
77 refGrad() = Zero;
78 valueFraction() = Zero;
79}
80
81
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_),
105)
106:
119)
120:
121 directionMixedFvPatchVectorField(pivpvf, iF),
122 phiName_(pivpvf.phiName_),
125{}
126
127
128// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129
131(
132 const fvPatchFieldMapper& m
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
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 auto& phip =
170 patch().lookupPatchField<surfaceScalarField>(phiName_);
171
172 valueFraction() += neg(phip)*(I - valueFraction());
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();
195}
196
197
198// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
199
200void Foam::fixedNormalInletOutletVelocityFvPatchVectorField::operator=
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);
209
210
211// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212
213namespace Foam
214{
216 (
218 fixedNormalInletOutletVelocityFvPatchVectorField
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...
@ MUST_READ
Reading required.
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
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.
fixedNormalInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
tmp< fvPatchVectorField > normalVelocity_
BC which provided the normal component of the velocity.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< fvPatchField< vector > > clone() const
Return a clone.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Switch fixTangentialInflow_
Set true to fix the tangential component for inflow.
const fvPatchVectorField & normalVelocity() const
Return the BC which provides the normal component of velocity.
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
Lookup type of boundary radiation properties.
Definition lookup.H:60
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,...
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< vector > fvPatchVectorField
dictionary dict
Foam::surfaceFields.