Loading...
Searching...
No Matches
fixedMultiPhaseHeatFluxFvPatchScalarField.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) 2015-2020 OpenFOAM Foundation
9 Copyright (C) 2020-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
33#include "phaseSystem.H"
35#include "ThermalDiffusivity.H"
37
38// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39
42(
43 const fvPatch& p,
45)
46:
47 fixedValueFvPatchScalarField(p, iF),
48 q_(p.size(), Zero),
49 relax_(1),
50 Tmin_(273)
51{}
52
53
56(
57 const fvPatch& p,
59 const dictionary& dict
60)
61:
62 fixedValueFvPatchScalarField(p, iF, dict),
63 q_("q", dict, p.size()),
64 relax_(dict.getOrDefault<scalar>("relax", 1)),
65 Tmin_(dict.getOrDefault<scalar>("Tmin", 273))
66{}
67
68
71(
73 const fvPatch& p,
75 const fvPatchFieldMapper& mapper
76)
77:
78 fixedValueFvPatchScalarField(psf, p, iF, mapper),
79 q_(psf.q_, mapper),
80 relax_(psf.relax_),
81 Tmin_(psf.Tmin_)
82{}
83
84
87(
89)
90:
91 fixedValueFvPatchScalarField(psf),
92 q_(psf.q_),
93 relax_(psf.relax_),
94 Tmin_(psf.Tmin_)
95{}
96
97
100(
103)
104:
105 fixedValueFvPatchScalarField(psf, iF),
106 q_(psf.q_),
107 relax_(psf.relax_),
108 Tmin_(psf.Tmin_)
109{}
110
111
112
113// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114
116{
117 if (updated())
118 {
119 return;
120 }
121
122 // Lookup the fluid model
123 const phaseSystem& fluid =
124 db().lookupObject<phaseSystem>("phaseProperties");
125
126 const scalarField& Tp = *this;
127
128 scalarField A(Tp.size(), Zero);
129 scalarField B(Tp.size(), Zero);
130 scalarField Q(Tp.size(), Zero);
131
132 forAll(fluid.phases(), phasei)
133 {
134 const phaseModel& phase = fluid.phases()[phasei];
135 const fluidThermo& thermo = phase.thermo();
136
138 phase.boundaryField()[patch().index()];
139
140 const fvPatchScalarField& T =
141 thermo.T().boundaryField()[patch().index()];
142
143 const scalarField kappaEff(phase.kappaEff(patch().index()));
144
145 if (debug)
146 {
147 const scalarField q0(T.snGrad()*alpha*kappaEff);
148 Q += q0;
149
150 auto limits = gMinMax(q0);
151
152 Info<< patch().name() << " " << phase.name()
153 << ": Heat flux "
154 << limits.min() << " - " << limits.max() << endl;
155 }
156
157 A += T.patchInternalField()*alpha*kappaEff*patch().deltaCoeffs();
158 B += alpha*kappaEff*patch().deltaCoeffs();
159 }
160
161 if (debug)
162 {
163 auto limits = gMinMax(Q);
164
165 Info<< patch().name() << " " << ": overall heat flux "
166 << limits.min() << " - " << limits.max() << " W/m2, power: "
167 << gWeightedSum(patch().magSf(), Q) << " W" << endl;
168 }
170 operator==((scalar(1) - relax_)*Tp + relax_*max(Tmin_,(q_ + A)/(B)));
171
172 fixedValueFvPatchScalarField::updateCoeffs();
173}
174
175
177(
178 const fvPatchFieldMapper& m
180{
181 fixedValueFvPatchScalarField::autoMap(m);
182 m(q_);
183}
184
185
187(
188 const fvPatchScalarField& ptf,
189 const labelList& addr
190)
191{
192 fixedValueFvPatchScalarField::rmap(ptf, addr);
193
196
197 q_.rmap(mptf.q_, addr);
198}
199
200
202{
204 os.writeEntry("relax", relax_);
205 q_.writeEntry("q", os);
208
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212namespace Foam
213{
215 (
217 fixedMultiPhaseHeatFluxFvPatchScalarField
218 );
219}
220
221
222// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Macros for easy insertion into run-time selection tables.
twoPhaseSystem & fluid
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Calculates a wall temperature that produces the specified overall wall heat flux across all the phase...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
fixedMultiPhaseHeatFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Fundamental fluid thermodynamic properties.
Definition fluidThermo.H:52
A FieldMapper for finite-volume patch fields.
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
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
Class to represent a system of phases and model interfacial transfers between them.
Definition phaseSystem.H:72
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const word & name() const
Definition phase.H:114
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
volScalarField & p
const volScalarField & T
auto limits
Definition setRDeltaT.H:186
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
label phasei
Definition pEqn.H:27
kappaEff
Definition TEqn.H:10
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fvPatchField< scalar > fvPatchScalarField
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299