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 | Copyright (C) 2012-2013 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 
27 #include "volFields.H"
29 #include "fvPatchFieldMapper.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
42  flowRate_(0),
43  alphaName_("none")
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
57  flowRate_(ptf.flowRate_),
58  alphaName_(ptf.alphaName_)
59 {}
60 
61 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
71  flowRate_(readScalar(dict.lookup("flowRate"))),
72  alphaName_(dict.lookup("alpha"))
73 {}
74 
75 
78 (
80 )
81 :
83  flowRate_(ptf.flowRate_),
84  alphaName_(ptf.alphaName_)
85 {}
86 
87 
90 (
93 )
94 :
96  flowRate_(ptf.flowRate_),
97  alphaName_(ptf.alphaName_)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
105 {
106  if (updated())
107  {
108  return;
109  }
110 
111  scalarField alphap =
112  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
113 
114  alphap = max(alphap, scalar(0));
115  alphap = min(alphap, scalar(1));
116 
117  // a simpler way of doing this would be nice
118  scalar avgU = -flowRate_/gSum(patch().magSf()*alphap);
119 
120  vectorField n(patch().nf());
121 
122  operator==(n*avgU*alphap);
123 
125 }
126 
127 
129 (
130  Ostream& os
131 ) const
132 {
134 
135  os.writeKeyword("flowRate") << flowRate_
136  << token::END_STATEMENT << nl;
137  os.writeKeyword("alpha") << alphaName_
138  << token::END_STATEMENT << nl;
139  writeEntry("value", os);
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 namespace Foam
146 {
148  (
151  );
152 }
153 
154 
155 // ************************************************************************* //
Foam::surfaceFields.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
#define readScalar
Definition: doubleScalar.C:38
Foam::fvPatchFieldMapper.
variableHeightFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
runTime write()
Type gSum(const FieldField< Field, Type > &f)
label n
dictionary dict
static const char nl
Definition: Ostream.H:260
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
volScalarField & p
Definition: createFields.H:51
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
This boundary condition provides a velocity boundary condition for multphase flow based on a user-spe...
Macros for easy insertion into run-time selection tables.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)