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-2023 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  refValue() = normalVelocity();
52  refGrad() = Zero;
53  valueFraction() = Zero;
54 }
55 
56 
59 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
66  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
67  phiName_(ptf.phiName_),
68  fixTangentialInflow_(ptf.fixTangentialInflow_),
69  normalVelocity_
70  (
71  fvPatchVectorField::New(ptf.normalVelocity(), p, iF, mapper)
72  )
73 {}
74 
75 
78 (
81 )
82 :
83  directionMixedFvPatchVectorField(pivpvf, iF),
84  phiName_(pivpvf.phiName_),
85  fixTangentialInflow_(pivpvf.fixTangentialInflow_),
86  normalVelocity_(pivpvf.normalVelocity().clone(iF))
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 (
94  const fvPatchVectorField& ptf,
95  const fvPatchFieldMapper& mapper
96 )
97 {
98  directionMixedFvPatchVectorField::map(ptf, mapper);
99 
101  refCast<const fixedNormalInletOutletVelocityFvPatchVectorField>(ptf);
102 
103  mapper(normalVelocity_.ref(), fniovptf.normalVelocity());
104 }
105 
106 
108 (
109  const fvPatchVectorField& ptf
110 )
111 {
112  directionMixedFvPatchVectorField::reset(ptf);
113 
115  refCast<const fixedNormalInletOutletVelocityFvPatchVectorField>(ptf);
116 
117  normalVelocity_->reset(fniovptf.normalVelocity());
118 }
119 
120 
122 {
123  if (updated())
124  {
125  return;
126  }
127 
128  normalVelocity_->evaluate();
129  refValue() = normalVelocity();
130 
131  valueFraction() = sqr(patch().nf());
132 
133  if (fixTangentialInflow_)
134  {
135  const fvsPatchField<scalar>& phip =
136  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
137 
138  valueFraction() += neg(phip)*(I - valueFraction());
139  }
140 
141  directionMixedFvPatchVectorField::updateCoeffs();
143 }
144 
145 
147 (
148  Ostream& os
149 )
150 const
151 {
153  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
154  writeEntry(os, "fixTangentialInflow", fixTangentialInflow_);
155  writeKeyword(os, "normalVelocity")
156  << nl << indent << token::BEGIN_BLOCK << nl << incrIndent;
157  normalVelocity_->write(os);
158  os << decrIndent << indent << token::END_BLOCK << endl;
159  writeEntry(os, "value", *this);
160 }
161 
162 
163 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
164 
165 void Foam::fixedNormalInletOutletVelocityFvPatchVectorField::operator=
166 (
167  const fvPatchField<vector>& pvf
168 )
169 {
170  tmp<vectorField> normalValue = transform(valueFraction(), refValue());
171  tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
172  fvPatchField<vector>::operator=(normalValue + transformGradValue);
173 }
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 namespace Foam
179 {
181  (
184  );
185 }
186 
187 // ************************************************************************* //
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: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:160
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 map(const fvPatchVectorField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
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.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
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:82
A class for managing temporary objects.
Definition: tmp.H:55
@ BEGIN_BLOCK
Definition: token.H:110
@ END_BLOCK
Definition: token.H:111
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:235
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
SurfaceField< scalar > surfaceScalarField
static const Identity< scalar > I
Definition: Identity.H:93
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
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:221
static const char nl
Definition: Ostream.H:260
dictionary dict
volScalarField & p
Foam::surfaceFields.