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-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 
27 #include "fieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  mixedFvPatchScalarField(p, iF, dict, false),
43  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
44  lowerBound_(dict.lookup<scalar>("lowerBound", dimless)),
45  upperBound_(dict.lookup<scalar>("upperBound", dimless))
46 {
47  this->refValue() = 0.0;
48 
49  if (dict.found("value"))
50  {
52  (
53  scalarField("value", iF.dimensions(), dict, p.size())
54  );
55  }
56  else
57  {
58  fvPatchScalarField::operator=(this->patchInternalField());
59  }
60 
61  this->refGrad() = 0.0;
62  this->valueFraction() = 0.0;
63 }
64 
65 
68 (
70  const fvPatch& p,
72  const fieldMapper& mapper
73 )
74 :
75  mixedFvPatchScalarField(ptf, p, iF, mapper),
76  phiName_(ptf.phiName_),
77  lowerBound_(ptf.lowerBound_),
78  upperBound_(ptf.upperBound_)
79 {}
80 
81 
84 (
87 )
88 :
89  mixedFvPatchScalarField(ptf, iF),
90  phiName_(ptf.phiName_),
91  lowerBound_(ptf.lowerBound_),
92  upperBound_(ptf.upperBound_)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
100  if (this->updated())
101  {
102  return;
103  }
104 
105  const fvsPatchField<scalar>& phip =
106  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
107 
108  scalarField alphap(this->patchInternalField());
109 
110 
111  forAll(phip, i)
112  {
113  if (phip[i] < -small)
114  {
115  if (alphap[i] < lowerBound_)
116  {
117  this->refValue()[i] = 0.0;
118  }
119  else if (alphap[i] > upperBound_)
120  {
121  this->refValue()[i] = 1.0;
122  }
123  else
124  {
125  this->refValue()[i] = alphap[i];
126  }
127 
128  this->valueFraction()[i] = 1.0;
129  }
130  else
131  {
132  this->refValue()[i] = 0.0;
133  this->valueFraction()[i] = 0.0;
134  }
135  }
136 
137  mixedFvPatchScalarField::updateCoeffs();
138 }
139 
140 
142 {
144  if (phiName_ != "phi")
145  {
146  writeEntry(os, "phi", phiName_);
147  }
148  writeEntry(os, "lowerBound", lowerBound_);
149  writeEntry(os, "upperBound", upperBound_);
150  writeEntry(os, "value", *this);
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 namespace Foam
157 {
159  (
162  );
163 }
164 
165 
166 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
const dimensionSet & dimensions() const
Return dimensions.
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
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
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
This boundary condition provides a phase fraction condition based on the local flow conditions,...
variableHeightFlowRateFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.