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 :
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 {
84  patchType() = dict.getOrDefault<word>("patchType", word::null);
85  this->refValue() = 0.0;
86 
87  if (dict.found("value"))
88  {
89  fvPatchScalarField::operator=
90  (
91  scalarField("value", dict, p.size())
92  );
93  }
94  else
95  {
96  fvPatchScalarField::operator=(this->patchInternalField());
97  }
98 
99  this->refGrad() = 0.0;
100  this->valueFraction() = 0.0;
101 }
102 
103 
106 (
108 )
109 :
110  mixedFvPatchScalarField(ptf),
111  phiName_(ptf.phiName_),
112  lowerBound_(ptf.lowerBound_),
113  upperBound_(ptf.upperBound_)
114 {}
115 
116 
119 (
122 )
123 :
124  mixedFvPatchScalarField(ptf, iF),
125  phiName_(ptf.phiName_),
126  lowerBound_(ptf.lowerBound_),
127  upperBound_(ptf.upperBound_)
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
134 {
135  if (this->updated())
136  {
137  return;
138  }
139 
140  const fvsPatchField<scalar>& phip =
141  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
142 
143  scalarField alphap(this->patchInternalField());
144 
145 
146  forAll(phip, i)
147  {
148  if (phip[i] < -SMALL)
149  {
150  if (alphap[i] < lowerBound_)
151  {
152  this->refValue()[i] = 0.0;
153  }
154  else if (alphap[i] > upperBound_)
155  {
156  this->refValue()[i] = 1.0;
157  }
158  else
159  {
160  this->refValue()[i] = alphap[i];
161  }
162 
163  this->valueFraction()[i] = 1.0;
164  }
165  else
166  {
167  this->refValue()[i] = 0.0;
168  this->valueFraction()[i] = 0.0;
169  }
170  }
171 
172  mixedFvPatchScalarField::updateCoeffs();
173 }
174 
175 
177 {
179  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
180  os.writeEntry("lowerBound", lowerBound_);
181  os.writeEntry("upperBound", upperBound_);
182  this->writeEntry("value", os);
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 namespace Foam
189 {
191  (
194  );
195 }
196 
197 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::variableHeightFlowRateFvPatchScalarField::lowerBound_
scalar lowerBound_
Lower bound for alpha1.
Definition: variableHeightFlowRateFvPatchField.H:116
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
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:248
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::variableHeightFlowRateFvPatchScalarField::upperBound_
scalar upperBound_
Upper bound for alpha1.
Definition: variableHeightFlowRateFvPatchField.H:119
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::variableHeightFlowRateFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: variableHeightFlowRateFvPatchField.C:133
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::variableHeightFlowRateFvPatchScalarField::phiName_
word phiName_
Name of flux field.
Definition: variableHeightFlowRateFvPatchField.H:113
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::variableHeightFlowRateFvPatchScalarField
This boundary condition provides a phase fraction condition based on the local flow conditions,...
Definition: variableHeightFlowRateFvPatchField.H:103
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::mixedFvPatchField< scalar >
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::variableHeightFlowRateFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: variableHeightFlowRateFvPatchField.C:176
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
variableHeightFlowRateFvPatchField.H
Foam::variableHeightFlowRateFvPatchScalarField::variableHeightFlowRateFvPatchScalarField
variableHeightFlowRateFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: variableHeightFlowRateFvPatchField.C:39
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54