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-------------------------------------------------------------------------------
11License
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
35Foam::variableHeightFlowRateInletVelocityFvPatchVectorField
36::variableHeightFlowRateInletVelocityFvPatchVectorField
37(
38 const fvPatch& p,
40)
41:
43 flowRate_(),
44 alphaName_("none")
45{}
46
47
48Foam::variableHeightFlowRateInletVelocityFvPatchVectorField
49::variableHeightFlowRateInletVelocityFvPatchVectorField
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
62Foam::variableHeightFlowRateInletVelocityFvPatchVectorField
63::variableHeightFlowRateInletVelocityFvPatchVectorField
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
77Foam::variableHeightFlowRateInletVelocityFvPatchVectorField
78::variableHeightFlowRateInletVelocityFvPatchVectorField
79(
81)
82:
84 flowRate_(ptf.flowRate_.clone()),
85 alphaName_(ptf.alphaName_)
86{}
87
88
89Foam::variableHeightFlowRateInletVelocityFvPatchVectorField
90::variableHeightFlowRateInletVelocityFvPatchVectorField
91(
94)
95:
97 flowRate_(ptf.flowRate_.clone()),
98 alphaName_(ptf.alphaName_)
99{}
100
101
102// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103
104void Foam::variableHeightFlowRateInletVelocityFvPatchVectorField
105::updateCoeffs()
106{
107 if (updated())
108 {
109 return;
110 }
111
112 scalarField alphap =
113 patch().lookupPatchField<volScalarField, scalar>(alphaName_);
114
115 alphap = max(alphap, scalar(0));
116 alphap = min(alphap, scalar(1));
117
118 const scalar t = db().time().timeOutputValue();
119 scalar flowRate = flowRate_->value(t);
120
121 // a simpler way of doing this would be nice
122 scalar avgU = -flowRate/gSum(patch().magSf()*alphap);
123
124 vectorField n(patch().nf());
125
126 operator==(n*avgU*alphap);
127
129}
130
131
133(
134 Ostream& os
135) const
136{
138 flowRate_->writeData(os);
139 os.writeEntry("alpha", alphaName_);
140 writeEntry("value", os);
141}
142
143
144// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145
146namespace Foam
147{
149 (
152 );
153}
154
155
156// ************************************************************************* //
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...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:126
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Lookup type of boundary radiation properties.
Definition: lookup.H:66
This boundary condition provides a velocity boundary condition for multphase flow based on a user-spe...
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict