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-2022 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 )
40 :
41  directionMixedFvPatchVectorField(p, iF),
42  phiName_("phi"),
43  fixTangentialInflow_(true),
44  normalVelocity_
45  (
46  fvPatchVectorField::New("fixedValue", p, iF)
47  )
48 {
49  refValue() = Zero;
50  refGrad() = Zero;
51  valueFraction() = Zero;
52 }
53 
54 
57 (
58  const fvPatch& p,
60  const dictionary& dict
61 )
62 :
63  directionMixedFvPatchVectorField(p, iF),
64  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
65  fixTangentialInflow_(dict.lookup("fixTangentialInflow")),
66  normalVelocity_
67  (
68  fvPatchVectorField::New(p, iF, dict.subDict("normalVelocity"))
69  )
70 {
71  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
72  refValue() = normalVelocity();
73  refGrad() = Zero;
74  valueFraction() = Zero;
75 }
76 
77 
80 (
82  const fvPatch& p,
84  const fvPatchFieldMapper& mapper
85 )
86 :
87  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
88  phiName_(ptf.phiName_),
89  fixTangentialInflow_(ptf.fixTangentialInflow_),
90  normalVelocity_
91  (
92  fvPatchVectorField::New(ptf.normalVelocity(), p, iF, mapper)
93  )
94 {}
95 
96 
99 (
102 )
103 :
104  directionMixedFvPatchVectorField(pivpvf, iF),
105  phiName_(pivpvf.phiName_),
106  fixTangentialInflow_(pivpvf.fixTangentialInflow_),
107  normalVelocity_(pivpvf.normalVelocity().clone(iF))
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 (
115  const fvPatchFieldMapper& m
116 )
117 {
118  directionMixedFvPatchVectorField::autoMap(m);
119  normalVelocity_->autoMap(m);
120 }
121 
122 
124 (
125  const fvPatchVectorField& ptf,
126  const labelList& addr
127 )
128 {
129  directionMixedFvPatchVectorField::rmap(ptf, addr);
130 
132  refCast<const fixedNormalInletOutletVelocityFvPatchVectorField>(ptf);
133 
134  normalVelocity_->rmap(fniovptf.normalVelocity(), addr);
135 }
136 
137 
139 (
140  const fvPatchVectorField& ptf
141 )
142 {
143  directionMixedFvPatchVectorField::reset(ptf);
144 
146  refCast<const fixedNormalInletOutletVelocityFvPatchVectorField>(ptf);
147 
148  normalVelocity_->reset(fniovptf.normalVelocity());
149 }
150 
151 
153 {
154  if (updated())
155  {
156  return;
157  }
158 
159  normalVelocity_->evaluate();
160  refValue() = normalVelocity();
161 
162  valueFraction() = sqr(patch().nf());
163 
164  if (fixTangentialInflow_)
165  {
166  const fvsPatchField<scalar>& phip =
167  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
168 
169  valueFraction() += neg(phip)*(I - valueFraction());
170  }
171 
172  directionMixedFvPatchVectorField::updateCoeffs();
174 }
175 
176 
178 (
179  Ostream& os
180 )
181 const
182 {
184  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
185  writeEntry(os, "fixTangentialInflow", fixTangentialInflow_);
186  writeKeyword(os, "normalVelocity")
187  << nl << indent << token::BEGIN_BLOCK << nl << incrIndent;
188  normalVelocity_->write(os);
189  os << decrIndent << indent << token::END_BLOCK << endl;
190  writeEntry(os, "value", *this);
191 }
192 
193 
194 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
195 
196 void Foam::fixedNormalInletOutletVelocityFvPatchVectorField::operator=
197 (
198  const fvPatchField<vector>& pvf
199 )
200 {
201  tmp<vectorField> normalValue = transform(valueFraction(), refValue());
202  tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
203  fvPatchField<vector>::operator=(normalValue + transformGradValue);
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 namespace Foam
210 {
212  (
215  );
216 }
217 
218 // ************************************************************************* //
Foam::surfaceFields.
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const fvPatchVectorField & normalVelocity() const
Return the BC which provides the normal component of velocity.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
This velocity inlet/outlet boundary condition combines a fixed normal component obtained from the "no...
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
volVectorField vectorField(fieldObject, mesh)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
static const Identity< scalar > I
Definition: Identity.H:93
A class for handling words, derived from string.
Definition: word.H:59
fixedNormalInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:258
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
tmp< fvPatchField< Type > > clone() const
Disallow clone without setting internal field reference.
Definition: fvPatchField.H:199
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864