fixedNormalInletOutletVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  directionMixedFvPatchVectorField(p, iF),
43  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
44  fixTangentialInflow_(dict.lookup("fixTangentialInflow")),
45  normalVelocity_
46  (
47  fvPatchVectorField::New(p, iF, dict.subDict("normalVelocity"))
48  )
49 {
51  (
52  vectorField("value", iF.dimensions(), dict, p.size())
53  );
54  refValue() = normalVelocity();
55  refGrad() = Zero;
56  valueFraction() = Zero;
57 }
58 
59 
62 (
64  const fvPatch& p,
66  const fieldMapper& mapper
67 )
68 :
69  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
70  phiName_(ptf.phiName_),
71  fixTangentialInflow_(ptf.fixTangentialInflow_),
72  normalVelocity_
73  (
74  fvPatchVectorField::New(ptf.normalVelocity(), p, iF, mapper)
75  )
76 {}
77 
78 
81 (
84 )
85 :
86  directionMixedFvPatchVectorField(pivpvf, iF),
87  phiName_(pivpvf.phiName_),
88  fixTangentialInflow_(pivpvf.fixTangentialInflow_),
89  normalVelocity_(pivpvf.normalVelocity().clone(iF))
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 (
97  const fvPatchVectorField& ptf,
98  const fieldMapper& mapper
99 )
100 {
101  directionMixedFvPatchVectorField::map(ptf, mapper);
102 
104  refCast<const fixedNormalInletOutletVelocityFvPatchVectorField>(ptf);
105 
106  mapper(normalVelocity_.ref(), fniovptf.normalVelocity());
107 }
108 
109 
111 (
112  const fvPatchVectorField& ptf
113 )
114 {
115  directionMixedFvPatchVectorField::reset(ptf);
116 
118  refCast<const fixedNormalInletOutletVelocityFvPatchVectorField>(ptf);
119 
120  normalVelocity_->reset(fniovptf.normalVelocity());
121 }
122 
123 
125 {
126  if (updated())
127  {
128  return;
129  }
130 
131  normalVelocity_->evaluate();
132  refValue() = normalVelocity();
133 
134  valueFraction() = sqr(patch().nf());
135 
136  if (fixTangentialInflow_)
137  {
138  const fvsPatchField<scalar>& phip =
139  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
140 
141  valueFraction() += neg(phip)*(I - valueFraction());
142  }
143 
144  directionMixedFvPatchVectorField::updateCoeffs();
146 }
147 
148 
150 (
151  Ostream& os
152 )
153 const
154 {
156  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
157  writeEntry(os, "fixTangentialInflow", fixTangentialInflow_);
158  writeKeyword(os, "normalVelocity")
159  << nl << indent << token::BEGIN_BLOCK << nl << incrIndent;
160  normalVelocity_->write(os);
161  os << decrIndent << indent << token::END_BLOCK << endl;
162  writeEntry(os, "value", *this);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
167 
168 void Foam::fixedNormalInletOutletVelocityFvPatchVectorField::operator=
169 (
170  const fvPatchField<vector>& pvf
171 )
172 {
173  tmp<vectorField> normalValue = transform(valueFraction(), refValue());
174  tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
175  fvPatchField<vector>::operator=(normalValue + transformGradValue);
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
187  );
188 }
189 
190 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This velocity inlet/outlet boundary condition combines a fixed normal component obtained from the "no...
fixedNormalInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
const fvPatchVectorField & normalVelocity() const
Return the BC which provides the normal component of velocity.
virtual void map(const fvPatchVectorField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
A class for managing temporary objects.
Definition: tmp.H:55
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
A class for handling words, derived from string.
Definition: word.H:62
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
SurfaceField< scalar > surfaceScalarField
static const Identity< scalar > I
Definition: Identity.H:93
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
T clone(const T &t)
Definition: List.H:55
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p
Foam::surfaceFields.