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 | Copyright (C) 2012 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 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  mixedFvPatchScalarField(ptf, p, iF, mapper),
62  phiName_(ptf.phiName_),
63  lowerBound_(ptf.lowerBound_),
64  upperBound_(ptf.upperBound_)
65 {}
66 
67 
70 (
71  const fvPatch& p,
73  const dictionary& dict
74 )
75 :
76  mixedFvPatchScalarField(p, iF),
77  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
78  lowerBound_(readScalar(dict.lookup("lowerBound"))),
79  upperBound_(readScalar(dict.lookup("upperBound")))
80 {
81  this->refValue() = 0.0;
82 
83  if (dict.found("value"))
84  {
85  fvPatchScalarField::operator=
86  (
87  scalarField("value", dict, p.size())
88  );
89  }
90  else
91  {
92  fvPatchScalarField::operator=(this->patchInternalField());
93  }
94 
95  this->refGrad() = 0.0;
96  this->valueFraction() = 0.0;
97 }
98 
99 
102 (
104 )
105 :
106  mixedFvPatchScalarField(ptf),
107  phiName_(ptf.phiName_),
108  lowerBound_(ptf.lowerBound_),
109  upperBound_(ptf.upperBound_)
110 {}
111 
112 
115 (
118 )
119 :
120  mixedFvPatchScalarField(ptf, iF),
121  phiName_(ptf.phiName_),
122  lowerBound_(ptf.lowerBound_),
123  upperBound_(ptf.upperBound_)
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
130 {
131  if (this->updated())
132  {
133  return;
134  }
135 
136  const fvsPatchField<scalar>& phip =
137  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
138 
139  scalarField alphap(this->patchInternalField());
140 
141 
142  forAll(phip, i)
143  {
144  if (phip[i] < -SMALL)
145  {
146  if (alphap[i] < lowerBound_)
147  {
148  this->refValue()[i] = 0.0;
149  }
150  else if (alphap[i] > upperBound_)
151  {
152  this->refValue()[i] = 1.0;
153  }
154  else
155  {
156  this->refValue()[i] = alphap[i];
157  }
158 
159  this->valueFraction()[i] = 1.0;
160  }
161  else
162  {
163  this->refValue()[i] = 0.0;
164  this->valueFraction()[i] = 0.0;
165  }
166  }
167 
168  mixedFvPatchScalarField::updateCoeffs();
169 }
170 
171 
173 {
175  if (phiName_ != "phi")
176  {
177  os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
178  }
179  os.writeKeyword("lowerBound") << lowerBound_ << token::END_STATEMENT << nl;
180  os.writeKeyword("upperBound") << upperBound_ << token::END_STATEMENT << nl;
181  this->writeEntry("value", os);
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 namespace Foam
188 {
190  (
193  );
194 }
195 
196 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::variableHeightFlowRateFvPatchScalarField::lowerBound_
scalar lowerBound_
Lower bound for alpha1.
Definition: variableHeightFlowRateFvPatchField.H:115
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::variableHeightFlowRateFvPatchScalarField::upperBound_
scalar upperBound_
Upper bound for alpha1.
Definition: variableHeightFlowRateFvPatchField.H:118
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::variableHeightFlowRateFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: variableHeightFlowRateFvPatchField.C:129
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::variableHeightFlowRateFvPatchScalarField::phiName_
word phiName_
Name of flux field.
Definition: variableHeightFlowRateFvPatchField.H:112
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::variableHeightFlowRateFvPatchScalarField
This boundary condition provides a phase fraction condition based on the local flow conditions,...
Definition: variableHeightFlowRateFvPatchField.H:102
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::mixedFvPatchField< scalar >
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::variableHeightFlowRateFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: variableHeightFlowRateFvPatchField.C:172
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
variableHeightFlowRateFvPatchField.H
Foam::variableHeightFlowRateFvPatchScalarField::variableHeightFlowRateFvPatchScalarField
variableHeightFlowRateFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: variableHeightFlowRateFvPatchField.C:36
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51