outletMappedUniformInletHeatAdditionFvPatchField.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) 2016-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "basicThermo.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedValueFvPatchScalarField(p, iF),
44  outletPatchName_(),
45  phiName_("phi"),
46  Q_(0),
47  TMin_(0),
48  TMax_(5000)
49 {}
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
62  outletPatchName_(ptf.outletPatchName_),
63  phiName_(ptf.phiName_),
64  Q_(ptf.Q_),
65  TMin_(ptf.TMin_),
66  TMax_(ptf.TMax_)
67 {}
68 
69 
72 (
73  const fvPatch& p,
75  const dictionary& dict
76 )
77 :
78  fixedValueFvPatchScalarField(p, iF, dict),
79  outletPatchName_(dict.get<word>("outletPatch")),
80  phiName_(dict.getOrDefault<word>("phi", "phi")),
81  Q_(dict.get<scalar>("Q")),
82  TMin_(dict.getOrDefault<scalar>("TMin", 0)),
83  TMax_(dict.getOrDefault<scalar>("TMax", 5000))
84 {}
85 
86 
87 
90 (
92 )
93 :
94  fixedValueFvPatchScalarField(ptf),
95  outletPatchName_(ptf.outletPatchName_),
96  phiName_(ptf.phiName_),
97  Q_(ptf.Q_),
98  TMin_(ptf.TMin_),
99  TMax_(ptf.TMax_)
100 {}
101 
102 
103 
106 (
109 )
110 :
111  fixedValueFvPatchScalarField(ptf, iF),
112  outletPatchName_(ptf.outletPatchName_),
113  phiName_(ptf.phiName_),
114  Q_(ptf.Q_),
115  TMin_(ptf.TMin_),
116  TMax_(ptf.TMax_)
117 {}
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
122 
124 {
125  if (this->updated())
126  {
127  return;
128  }
129 
130  const volScalarField& vsf =
131  (
132  dynamic_cast<const volScalarField&>(this->internalField())
133  );
134 
135  const fvPatch& fvp = this->patch();
136 
137  label outletPatchID =
138  fvp.patch().boundaryMesh().findPatchID(outletPatchName_);
139 
140  if (outletPatchID < 0)
141  {
143  << "Unable to find outlet patch " << outletPatchName_
144  << abort(FatalError);
145  }
146 
147  const fvPatch& outletPatch = fvp.boundaryMesh()[outletPatchID];
148 
149  const fvPatchField<scalar>& outletPatchField =
150  vsf.boundaryField()[outletPatchID];
151 
152  const surfaceScalarField& phi =
153  db().lookupObject<surfaceScalarField>(phiName_);
154 
155  const scalarField& outletPatchPhi = phi.boundaryField()[outletPatchID];
156  scalar sumOutletPatchPhi = gSum(outletPatchPhi);
157 
158  if (sumOutletPatchPhi > SMALL)
159  {
160  const basicThermo& thermo =
161  db().lookupObject<basicThermo>(basicThermo::dictName);
162 
163  const scalarField& pp = thermo.p().boundaryField()[outletPatchID];
164  const scalarField& pT = thermo.T().boundaryField()[outletPatchID];
165 
166  scalar averageOutletField =
167  gSum(outletPatchPhi*outletPatchField)/sumOutletPatchPhi;
168 
169  const scalarField Cpf(thermo.Cp(pp, pT, outletPatchID));
170 
171  scalar totalPhiCp = gSum(outletPatchPhi)*gAverage(Cpf);
172 
173  operator==(min(max(averageOutletField + Q_/totalPhiCp, TMin_), TMax_));
174  }
175  else
176  {
177  scalar averageOutletField =
178  gSum(outletPatch.magSf()*outletPatchField)
179  /gSum(outletPatch.magSf());
180 
181  operator==(averageOutletField);
182  }
183 
184  fixedValueFvPatchScalarField::updateCoeffs();
185 }
186 
187 
189 (
190  Ostream& os
191 ) const
192 {
194  os.writeEntry("outletPatch", outletPatchName_);
195  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
196  os.writeEntry("Q", Q_);
197  os.writeEntry("TMin", TMin_);
198  os.writeEntry("TMax", TMax_);
199 
200  this->writeEntry("value", os);
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 namespace Foam
207 {
209  (
212  );
213 }
214 
215 
216 // ************************************************************************* //
Foam::fvPatchField< scalar >
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
p
volScalarField & p
Definition: createFieldRefs.H:8
basicThermo.H
Foam::outletMappedUniformInletHeatAdditionFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: outletMappedUniformInletHeatAdditionFvPatchField.C:189
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::outletMappedUniformInletHeatAdditionFvPatchField::outletMappedUniformInletHeatAdditionFvPatchField
outletMappedUniformInletHeatAdditionFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: outletMappedUniformInletHeatAdditionFvPatchField.C:38
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::fvPatch::boundaryMesh
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:203
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
Foam::outletMappedUniformInletHeatAdditionFvPatchField
This temperature boundary condition averages the temperature over the "outlet" patch specified by nam...
Definition: outletMappedUniformInletHeatAdditionFvPatchField.H:113
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:765
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::outletMappedUniformInletHeatAdditionFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: outletMappedUniformInletHeatAdditionFvPatchField.C:123
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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
outletMappedUniformInletHeatAdditionFvPatchField.H
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:161
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
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)
Foam::fvPatch::magSf
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:156
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60