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-2024 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  const dictionary& dict
38 )
39 :
41  flowRate_
42  (
43  Function1<scalar>::New
44  (
45  "flowRate",
46  db().time().userUnits(),
48  dict
49  )
50  ),
51  alphaName_(dict.lookup("alpha"))
52 {}
53 
54 
57 (
59  const fvPatch& p,
61  const fieldMapper& mapper
62 )
63 :
64  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
65  flowRate_(ptf.flowRate_, false),
66  alphaName_(ptf.alphaName_)
67 {}
68 
69 
72 (
75 )
76 :
78  flowRate_(ptf.flowRate_, false),
79  alphaName_(ptf.alphaName_)
80 {}
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
87 {
88  if (updated())
89  {
90  return;
91  }
92 
93  scalarField alphap =
94  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
95 
96  alphap = max(alphap, scalar(0));
97  alphap = min(alphap, scalar(1));
98 
99  scalar flowRate = flowRate_->value(db().time().value());
100 
101  // a simpler way of doing this would be nice
102  scalar avgU = -flowRate/gSum(patch().magSf()*alphap);
103 
104  vectorField n(patch().nf());
105 
106  operator==(n*avgU*alphap);
107 
109 }
110 
111 
113 (
114  Ostream& os
115 ) const
116 {
118  writeEntry(os, db().time().userUnits(), dimVolumetricFlux, flowRate_());
119  writeEntry(os, "alpha", alphaName_);
120  writeEntry(os, "value", *this);
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 namespace Foam
127 {
129  (
132  );
133 }
134 
135 
136 // ************************************************************************* //
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...
Run-time selectable general function of one variable.
Definition: Function1.H:125
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a velocity boundary condition for multiphase flow based on a user-sp...
variableHeightFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimVolumetricFlux
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p