variableHeightFlowRateInletVelocityFvPatchVectorField.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2021 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 
31 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
43  flowRate_(),
44  alphaName_("none")
45 {}
46 
47 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
57  flowRate_(Function1<scalar>::New("flowRate", dict, &db())),
58  alphaName_(dict.lookup("alpha"))
59 {}
60 
61 
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
72  flowRate_(ptf.flowRate_.clone()),
73  alphaName_(ptf.alphaName_)
74 {}
75 
76 
79 (
81 )
82 :
84  flowRate_(ptf.flowRate_.clone()),
85  alphaName_(ptf.alphaName_)
86 {}
87 
88 
91 (
94 )
95 :
97  flowRate_(ptf.flowRate_.clone()),
98  alphaName_(ptf.alphaName_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 {
106  if (updated())
107  {
108  return;
109  }
110 
111  const auto& mesh = patch().boundaryMesh().mesh();
112  auto& alpha = mesh.lookupObjectRef<volScalarField>(alphaName_);
113 
114  // Update alpha boundary (if needed) due to mesh changes
115  if (!mesh.upToDatePoints(alpha))
116  {
117  auto& alphabf = alpha.boundaryFieldRef();
118 
119  if (!alphabf[patch().index()].updated())
120  {
121  DebugInfo<< "Updating alpha BC due to mesh changes" << endl;
122  alphabf.evaluateSelected(labelList({ patch().index() }));
123  }
124  }
125 
126  scalarField alphap(alpha.boundaryField()[patch().index()]);
127 
128  alphap = max(alphap, scalar(0));
129  alphap = min(alphap, scalar(1));
130 
131  const scalar t = db().time().timeOutputValue();
132  scalar flowRate = flowRate_->value(t);
133 
134  // a simpler way of doing this would be nice
135  scalar avgU = -flowRate/gSum(patch().magSf()*alphap);
136 
137  vectorField n(patch().nf());
138 
139  operator==(n*avgU*alphap);
140 
142 }
143 
144 
146 (
147  Ostream& os
148 ) const
149 {
151  flowRate_->writeData(os);
152  os.writeEntry("alpha", alphaName_);
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 namespace Foam
160 {
162  (
164  variableHeightFlowRateInletVelocityFvPatchVectorField
165  );
166 }
167 
168 
169 // ************************************************************************* //
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
fvPatchField< vector > fvPatchVectorField
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
This boundary condition provides a velocity boundary condition for multphase flow based on a user-spe...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
variableHeightFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
#define DebugInfo
Report an information message using Foam::Info.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
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.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.