variableHeightFlowRateFvPatchField.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 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  mixedFvPatchField<scalar>(p, iF),
45  phiName_("phi"),
46  lowerBound_(0.0),
47  upperBound_(1.0)
48 {
49  this->refValue() = 0.0;
50  this->refGrad() = 0.0;
51  this->valueFraction() = 0.0;
52 }
53 
54 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  mixedFvPatchScalarField(ptf, p, iF, mapper),
65  phiName_(ptf.phiName_),
66  lowerBound_(ptf.lowerBound_),
67  upperBound_(ptf.upperBound_)
68 {}
69 
70 
73 (
74  const fvPatch& p,
76  const dictionary& dict
77 )
78 :
79  mixedFvPatchScalarField(p, iF),
80  phiName_(dict.getOrDefault<word>("phi", "phi")),
81  lowerBound_(dict.get<scalar>("lowerBound")),
82  upperBound_(dict.get<scalar>("upperBound"))
83 {
85 
86  if (!this->readValueEntry(dict))
87  {
88  // Fallback: set to the internal field
89  this->extrapolateInternal();
90  }
91 
92  this->refValue() = 0.0;
93  this->refGrad() = 0.0;
94  this->valueFraction() = 0.0;
95 }
96 
97 
100 (
102 )
103 :
104  mixedFvPatchScalarField(ptf),
105  phiName_(ptf.phiName_),
106  lowerBound_(ptf.lowerBound_),
107  upperBound_(ptf.upperBound_)
108 {}
109 
110 
113 (
116 )
117 :
118  mixedFvPatchScalarField(ptf, iF),
119  phiName_(ptf.phiName_),
120  lowerBound_(ptf.lowerBound_),
121  upperBound_(ptf.upperBound_)
122 {}
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
128 {
129  if (this->updated())
130  {
131  return;
132  }
133 
134  const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
135 
136  scalarField alphap(this->patchInternalField());
137 
138 
139  forAll(phip, i)
140  {
141  if (phip[i] < -SMALL)
142  {
143  if (alphap[i] < lowerBound_)
144  {
145  this->refValue()[i] = 0.0;
146  }
147  else if (alphap[i] > upperBound_)
148  {
149  this->refValue()[i] = 1.0;
150  }
151  else
152  {
153  this->refValue()[i] = alphap[i];
154  }
155 
156  this->valueFraction()[i] = 1.0;
157  }
158  else
159  {
160  this->refValue()[i] = 0.0;
161  this->valueFraction()[i] = 0.0;
162  }
163  }
164 
165  mixedFvPatchScalarField::updateCoeffs();
166 }
167 
168 
170 {
172  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
173  os.writeEntry("lowerBound", lowerBound_);
174  os.writeEntry("upperBound", upperBound_);
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
186  variableHeightFlowRateFvPatchScalarField
187  );
188 }
189 
190 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void readDict(const dictionary &dict)
Read dictionary entries.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
This boundary condition provides a phase fraction condition based on the local flow conditions...
A FieldMapper for finite-volume patch fields.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
OBJstream os(runTime.globalPath()/outputName)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
variableHeightFlowRateFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.