outletPhaseMeanVelocityFvPatchVectorField.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 #include "volFields.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
45  Umean_(0),
46  alphaName_("none")
47 {
48  refValue() = Zero;
49  refGrad() = Zero;
50  valueFraction() = 0.0;
51 }
52 
53 
56 (
58  const fvPatch& p,
60  const fvPatchFieldMapper& mapper
61 )
62 :
63  mixedFvPatchField<vector>(ptf, p, iF, mapper),
64  Umean_(ptf.Umean_),
65  alphaName_(ptf.alphaName_)
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
78  Umean_(dict.get<scalar>("Umean")),
79  alphaName_(dict.lookup("alpha"))
80 {
82 
83  refValue() = Zero;
84  refGrad() = Zero;
85  valueFraction() = 0;
86 
87  if (!this->readValueEntry(dict))
88  {
89  this->extrapolateInternal();
90  }
91 }
92 
93 
96 (
97  const outletPhaseMeanVelocityFvPatchVectorField& ptf
98 )
99 :
100  mixedFvPatchField<vector>(ptf),
101  Umean_(ptf.Umean_),
102  alphaName_(ptf.alphaName_)
103 {}
104 
105 
108 (
111 )
112 :
113  mixedFvPatchField<vector>(ptf, iF),
114  Umean_(ptf.Umean_),
115  alphaName_(ptf.alphaName_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
123  if (updated())
124  {
125  return;
126  }
127 
128  scalarField alphap
129  (
130  patch().lookupPatchField<volScalarField>(alphaName_)
131  );
132 
133  alphap = max(alphap, scalar(0));
134  alphap = min(alphap, scalar(1));
135 
136  // Get the patchInternalField (zero-gradient field)
137  vectorField Uzg(patchInternalField());
138 
139  // Calculate the phase mean zero-gradient velocity
140  scalar Uzgmean =
141  gSum(alphap*(patch().Sf() & Uzg))
142  /gSum(alphap*patch().magSf());
143 
144  // Set the refValue and valueFraction to adjust the boundary field
145  // such that the phase mean is Umean_
146  if (Uzgmean >= Umean_)
147  {
148  refValue() = Zero;
149  valueFraction() = 1.0 - Umean_/Uzgmean;
150  }
151  else
152  {
153  refValue() = (Umean_ + Uzgmean)*patch().nf();
154  valueFraction() = 1.0 - Uzgmean/Umean_;
155  }
156 
158 }
159 
160 
162 (
163  Ostream& os
164 ) const
165 {
167 
168  os.writeEntry("Umean", Umean_);
169  os.writeEntry("alpha", alphaName_);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 namespace Foam
177 {
179  (
181  outletPhaseMeanVelocityFvPatchVectorField
182  );
183 }
184 
185 
186 // ************************************************************************* //
outletPhaseMeanVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
dictionary dict
fvPatchField< vector > fvPatchVectorField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
This boundary condition adjusts the velocity for the given phase to achieve the specified mean thus c...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void readDict(const dictionary &dict)
Read dictionary entries.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
Type gSum(const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A FieldMapper for finite-volume patch fields.
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127