pressureInletOutletVelocityFvPatchVectorField.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) 2011-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 {
46  (
47  vectorField("value", iF.dimensions(), dict, p.size())
48  );
49 
50  if (dict.found("tangentialVelocity"))
51  {
52  tangentialVelocity_ =
54  (
55  "tangentialVelocity",
56  db().time().userUnits(),
58  dict
59  );
60  }
61 
62  refValue() = Zero;
63  refGrad() = Zero;
64  valueFraction() = Zero;
65 }
66 
67 
70 (
72  const fvPatch& p,
74  const fieldMapper& mapper
75 )
76 :
77  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
78  phiName_(ptf.phiName_),
79  tangentialVelocity_(ptf.tangentialVelocity_, false)
80 {}
81 
82 
85 (
88 )
89 :
90  directionMixedFvPatchVectorField(pivpvf, iF),
91  phiName_(pivpvf.phiName_),
92  tangentialVelocity_(pivpvf.tangentialVelocity_, false)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
100  if (updated())
101  {
102  return;
103  }
104 
105  if (tangentialVelocity_.valid())
106  {
107  const vector tangentialVelocity =
108  tangentialVelocity_->value(db().time().value());
109  const vectorField n(patch().nf());
110  refValue() = tangentialVelocity - n*(n & tangentialVelocity);
111  }
112 
113  const fvsPatchField<scalar>& phip =
114  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
115 
116  valueFraction() = neg(phip)*(I - sqr(patch().nf()));
117 
118  directionMixedFvPatchVectorField::updateCoeffs();
120 }
121 
122 
124 (
125  Ostream& os
126 )
127 const
128 {
130  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
131  if (tangentialVelocity_.valid())
132  {
133  writeEntry
134  (
135  os,
136  db().time().userUnits(),
137  dimVelocity,
138  tangentialVelocity_()
139  );
140  }
141  writeEntry(os, "value", *this);
142 }
143 
144 
145 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
146 
147 void Foam::pressureInletOutletVelocityFvPatchVectorField::operator=
148 (
149  const fvPatchField<vector>& pvf
150 )
151 {
152  tmp<vectorField> normalValue = transform(valueFraction(), refValue());
153  tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
154  fvPatchField<vector>::operator=(normalValue + transformGradValue);
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 namespace Foam
161 {
163  (
166  );
167 }
168 
169 // ************************************************************************* //
label n
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.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
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
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
Velocity inlet/outlet boundary condition for patches where the pressure is specified in some manner,...
pressureInletOutletVelocityFvPatchVectorField(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.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedSymmTensor sqr(const dimensionedVector &dv)
SurfaceField< scalar > surfaceScalarField
static const Identity< scalar > I
Definition: Identity.H:93
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimVelocity
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.