Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 Qptr_(nullptr),
45 outletPatchName_(),
46 phiName_("phi"),
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 Qptr_(ptf.Qptr_.clone()),
63 outletPatchName_(ptf.outletPatchName_),
64 phiName_(ptf.phiName_),
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 Qptr_(Function1<scalar>::New("Q", dict, &db())),
80 outletPatchName_(dict.get<word>("outletPatch")),
81 phiName_(dict.getOrDefault<word>("phi", "phi")),
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 Qptr_(ptf.Qptr_.clone()),
96 outletPatchName_(ptf.outletPatchName_),
97 phiName_(ptf.phiName_),
98 TMin_(ptf.TMin_),
99 TMax_(ptf.TMax_)
100{}
101
102
103
106(
109)
110:
111 fixedValueFvPatchScalarField(ptf, iF),
112 Qptr_(ptf.Qptr_.clone()),
113 outletPatchName_(ptf.outletPatchName_),
114 phiName_(ptf.phiName_),
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 const 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 gWeightedSum(outletPatchPhi, outletPatchField)/sumOutletPatchPhi;
168
169 // Calculate Q as a function of average outlet temperature
170 const scalar Q = Qptr_->value(averageOutletField);
171
172 const scalarField Cpf(thermo.Cp(pp, pT, outletPatchID));
173
174 scalar totalPhiCp = sumOutletPatchPhi*gAverage(Cpf);
175
176 operator==(clamp(averageOutletField + Q/totalPhiCp, TMin_, TMax_));
177 }
178 else
179 {
180 scalar averageOutletField =
181 gWeightedAverage(outletPatch.magSf(), outletPatchField);
182
183 operator==(averageOutletField);
184 }
185
186 fixedValueFvPatchScalarField::updateCoeffs();
187}
188
189
191(
192 Ostream& os
193) const
194{
196 os.writeEntry("outletPatch", outletPatchName_);
197 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
198
199 Qptr_->writeData(os);
200
201 os.writeEntry("TMin", TMin_);
202 os.writeEntry("TMax", TMax_);
203
206
207
208// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209
210namespace Foam
211{
213 (
215 outletMappedUniformInletHeatAdditionFvPatchField
216 );
217}
218
219
220// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors).
Definition Field.C:745
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition fvPatch.H:268
const scalarField & magSf() const
Return face area magnitudes, like the fvMesh::magSf() method.
Definition fvPatch.C:131
This temperature boundary condition averages the temperature over the "outlet" patch specified by nam...
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
outletMappedUniformInletHeatAdditionFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition polyPatch.C:295
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
fvPatchField< scalar > fvPatchScalarField
dictionary dict
psiReactionThermo & thermo
Foam::surfaceFields.