variableHeightFlowRateInletVelocityFvPatchVectorField.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) 2012-2021 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 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
34 (
35  const fvPatch& p,
37 )
38 :
40  flowRate_(),
41  alphaName_("none")
42 {}
43 
44 
47 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
54  flowRate_(Function1<scalar>::New("flowRate", dict)),
55  alphaName_(dict.lookup("alpha"))
56 {}
57 
58 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
69  flowRate_(ptf.flowRate_, false),
70  alphaName_(ptf.alphaName_)
71 {}
72 
73 
76 (
79 )
80 :
82  flowRate_(ptf.flowRate_, false),
83  alphaName_(ptf.alphaName_)
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
91 {
92  if (updated())
93  {
94  return;
95  }
96 
97  scalarField alphap =
98  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
99 
100  alphap = max(alphap, scalar(0));
101  alphap = min(alphap, scalar(1));
102 
103  const scalar t = db().time().userTimeValue();
104  scalar flowRate = flowRate_->value(t);
105 
106  // a simpler way of doing this would be nice
107  scalar avgU = -flowRate/gSum(patch().magSf()*alphap);
108 
109  vectorField n(patch().nf());
110 
111  operator==(n*avgU*alphap);
112 
114 }
115 
116 
118 (
119  Ostream& os
120 ) const
121 {
123  writeEntry(os, flowRate_());
124  writeEntry(os, "alpha", alphaName_);
125  writeEntry(os, "value", *this);
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 
131 namespace Foam
132 {
134  (
137  );
138 }
139 
140 
141 // ************************************************************************* //
Run-time selectable general function of one variable.
Definition: Function1.H:52
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
This boundary condition provides a velocity boundary condition for multiphase flow based on a user-sp...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::fvPatchFieldMapper.
variableHeightFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
label n
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864