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-2023 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_(Function1<scalar>::New("flowRate", dict)),
42  alphaName_(dict.lookup("alpha"))
43 {}
44 
45 
48 (
50  const fvPatch& p,
52  const fvPatchFieldMapper& mapper
53 )
54 :
55  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
56  flowRate_(ptf.flowRate_, false),
57  alphaName_(ptf.alphaName_)
58 {}
59 
60 
63 (
66 )
67 :
69  flowRate_(ptf.flowRate_, false),
70  alphaName_(ptf.alphaName_)
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
78 {
79  if (updated())
80  {
81  return;
82  }
83 
84  scalarField alphap =
85  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
86 
87  alphap = max(alphap, scalar(0));
88  alphap = min(alphap, scalar(1));
89 
90  const scalar t = db().time().userTimeValue();
91  scalar flowRate = flowRate_->value(t);
92 
93  // a simpler way of doing this would be nice
94  scalar avgU = -flowRate/gSum(patch().magSf()*alphap);
95 
96  vectorField n(patch().nf());
97 
98  operator==(n*avgU*alphap);
99 
101 }
102 
103 
105 (
106  Ostream& os
107 ) const
108 {
110  writeEntry(os, flowRate_());
111  writeEntry(os, "alpha", alphaName_);
112  writeEntry(os, "value", *this);
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 namespace Foam
119 {
121  (
124  );
125 }
126 
127 
128 // ************************************************************************* //
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:64
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:160
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
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 > &)
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:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p