variableHeightFlowRateFvPatchField.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 
27 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
42  phiName_("phi"),
43  lowerBound_(0.0),
44  upperBound_(1.0)
45 {
46  this->refValue() = 0.0;
47  this->refGrad() = 0.0;
48  this->valueFraction() = 0.0;
49 }
50 
51 
54 (
55  const fvPatch& p,
57  const dictionary& dict
58 )
59 :
60  mixedFvPatchScalarField(p, iF),
61  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
62  lowerBound_(dict.lookup<scalar>("lowerBound")),
63  upperBound_(dict.lookup<scalar>("upperBound"))
64 {
65  this->refValue() = 0.0;
66 
67  if (dict.found("value"))
68  {
69  fvPatchScalarField::operator=
70  (
71  scalarField("value", dict, p.size())
72  );
73  }
74  else
75  {
76  fvPatchScalarField::operator=(this->patchInternalField());
77  }
78 
79  this->refGrad() = 0.0;
80  this->valueFraction() = 0.0;
81 }
82 
83 
86 (
88  const fvPatch& p,
90  const fvPatchFieldMapper& mapper
91 )
92 :
93  mixedFvPatchScalarField(ptf, p, iF, mapper),
94  phiName_(ptf.phiName_),
95  lowerBound_(ptf.lowerBound_),
96  upperBound_(ptf.upperBound_)
97 {}
98 
99 
102 (
105 )
106 :
107  mixedFvPatchScalarField(ptf, iF),
108  phiName_(ptf.phiName_),
109  lowerBound_(ptf.lowerBound_),
110  upperBound_(ptf.upperBound_)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  if (this->updated())
119  {
120  return;
121  }
122 
123  const fvsPatchField<scalar>& phip =
124  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
125 
126  scalarField alphap(this->patchInternalField());
127 
128 
129  forAll(phip, i)
130  {
131  if (phip[i] < -small)
132  {
133  if (alphap[i] < lowerBound_)
134  {
135  this->refValue()[i] = 0.0;
136  }
137  else if (alphap[i] > upperBound_)
138  {
139  this->refValue()[i] = 1.0;
140  }
141  else
142  {
143  this->refValue()[i] = alphap[i];
144  }
145 
146  this->valueFraction()[i] = 1.0;
147  }
148  else
149  {
150  this->refValue()[i] = 0.0;
151  this->valueFraction()[i] = 0.0;
152  }
153  }
154 
155  mixedFvPatchScalarField::updateCoeffs();
156 }
157 
158 
160 {
162  if (phiName_ != "phi")
163  {
164  writeEntry(os, "phi", phiName_);
165  }
166  writeEntry(os, "lowerBound", lowerBound_);
167  writeEntry(os, "upperBound", upperBound_);
168  writeEntry(os, "value", *this);
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 namespace Foam
175 {
177  (
180  );
181 }
182 
183 
184 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:62
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:260
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from string.
Definition: word.H:59
This boundary condition provides a phase fraction condition based on the local flow conditions...
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:156
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
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.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
variableHeightFlowRateFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844