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"
32 #include "fvPatchFieldMapper.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 {
81  patchType() = dict.getOrDefault<word>("patchType", word::null);
82 
83  refValue() = Zero;
84  refGrad() = Zero;
85  valueFraction() = 0.0;
86 
87  if (dict.found("value"))
88  {
89  fvPatchVectorField::operator=
90  (
91  vectorField("value", dict, p.size())
92  );
93  }
94  else
95  {
96  fvPatchVectorField::operator=(patchInternalField());
97  }
98 }
99 
100 
103 (
105 )
106 :
108  Umean_(ptf.Umean_),
109  alphaName_(ptf.alphaName_)
110 {}
111 
112 
115 (
118 )
119 :
120  mixedFvPatchField<vector>(ptf, iF),
121  Umean_(ptf.Umean_),
122  alphaName_(ptf.alphaName_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 {
130  if (updated())
131  {
132  return;
133  }
134 
135  scalarField alphap =
136  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
137 
138  alphap = max(alphap, scalar(0));
139  alphap = min(alphap, scalar(1));
140 
141  // Get the patchInternalField (zero-gradient field)
142  vectorField Uzg(patchInternalField());
143 
144  // Calculate the phase mean zero-gradient velocity
145  scalar Uzgmean =
146  gSum(alphap*(patch().Sf() & Uzg))
147  /gSum(alphap*patch().magSf());
148 
149  // Set the refValue and valueFraction to adjust the boundary field
150  // such that the phase mean is Umean_
151  if (Uzgmean >= Umean_)
152  {
153  refValue() = Zero;
154  valueFraction() = 1.0 - Umean_/Uzgmean;
155  }
156  else
157  {
158  refValue() = (Umean_ + Uzgmean)*patch().nf();
159  valueFraction() = 1.0 - Uzgmean/Umean_;
160  }
161 
163 }
164 
165 
167 (
168  Ostream& os
169 ) const
170 {
172 
173  os.writeEntry("Umean", Umean_);
174  os.writeEntry("alpha", alphaName_);
175  writeEntry("value", os);
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
187  );
188 }
189 
190 
191 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
outletPhaseMeanVelocityFvPatchVectorField.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::outletPhaseMeanVelocityFvPatchVectorField::outletPhaseMeanVelocityFvPatchVectorField
outletPhaseMeanVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: outletPhaseMeanVelocityFvPatchVectorField.C:39
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::mixedFvPatchField< vector >
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::outletPhaseMeanVelocityFvPatchVectorField
This boundary condition adjusts the velocity for the given phase to achieve the specified mean thus c...
Definition: outletPhaseMeanVelocityFvPatchVectorField.H:94
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::outletPhaseMeanVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: outletPhaseMeanVelocityFvPatchVectorField.C:128
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::outletPhaseMeanVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: outletPhaseMeanVelocityFvPatchVectorField.C:167
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54