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-2019 OpenFOAM Foundation
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 
29 #include "fvPatchFieldMapper.H"
31 
32 #include "phaseSystem.H"
34 #include "ThermalDiffusivity.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44 )
45 :
46  fixedValueFvPatchScalarField(p, iF),
47  q_(p.size(), Zero),
48  relax_(1.0),
49  Tmin_(0.0)
50 {}
51 
52 
55 (
56  const fvPatch& p,
58  const dictionary& dict
59 )
60 :
61  fixedValueFvPatchScalarField(p, iF, dict),
62  q_("q", dict, p.size()),
63  relax_(dict.lookupOrDefault<scalar>("relax", 1.0)),
64  Tmin_(dict.lookupOrDefault<scalar>("Tmin", 273))
65 {}
66 
67 
70 (
72  const fvPatch& p,
74  const fvPatchFieldMapper& mapper
75 )
76 :
77  fixedValueFvPatchScalarField(psf, p, iF, mapper),
78  q_(psf.q_, mapper),
79  relax_(psf.relax_),
80  Tmin_(psf.Tmin_)
81 {}
82 
83 
86 (
88 )
89 :
90  fixedValueFvPatchScalarField(psf),
91  q_(psf.q_),
92  relax_(psf.relax_),
93  Tmin_(psf.Tmin_)
94 {}
95 
96 
99 (
102 )
103 :
104  fixedValueFvPatchScalarField(psf, iF),
105  q_(psf.q_),
106  relax_(psf.relax_),
107  Tmin_(psf.Tmin_)
108 {}
109 
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  if (updated())
117  {
118  return;
119  }
120 
121  // Lookup the fluid model
122  const phaseSystem& fluid =
123  (
124  db().lookupObject<phaseSystem>("phaseProperties")
125  );
126 
127  const scalarField& Tp = *this;
128 
129  scalarField A(Tp.size(), Zero);
130  scalarField B(Tp.size(), Zero);
131  scalarField Q(Tp.size(), Zero);
132 
133  forAll(fluid.phases(), phasei)
134  {
135  const phaseModel& phase = fluid.phases()[phasei];
136  const fluidThermo& thermo = phase.thermo();
137 
138  const fvPatchScalarField& alpha =
139  phase.boundaryField()[patch().index()];
140 
141  const fvPatchScalarField& T =
142  thermo.T().boundaryField()[patch().index()];
143 
144  const scalarField kappaEff(phase.kappaEff(patch().index()));
145 
146  if (debug)
147  {
148  scalarField q0(T.snGrad()*alpha*kappaEff);
149  Q += q0;
150 
151  Info<< patch().name() << " " << phase.name()
152  << ": Heat flux " << gMin(q0) << " - " << gMax(q0) << endl;
153  }
154 
155  A += T.patchInternalField()*alpha*kappaEff*patch().deltaCoeffs();
156  B += alpha*kappaEff*patch().deltaCoeffs();
157  }
158 
159  if (debug)
160  {
161  Info<< patch().name() << " " << ": overall heat flux "
162  << gMin(Q) << " - " << gMax(Q) << " W/m2, power: "
163  << gSum(patch().magSf()*Q) << " W" << endl;
164  }
165 
166  operator==((1 - relax_)*Tp + relax_*max(Tmin_,(q_ + A)/(B)));
167 
168  fixedValueFvPatchScalarField::updateCoeffs();
169 }
170 
171 
173 {
175  os.writeEntry("relax", relax_);
176  q_.writeEntry("q", os);
177  writeEntry("value", os);
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 namespace Foam
184 {
186  (
189  );
190 }
191 
192 
193 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:57
p
volScalarField & p
Definition: createFieldRefs.H:8
compressibleTurbulenceModel.H
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
ThermalDiffusivity.H
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
kappaEff
kappaEff
Definition: TEqn.H:7
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:52
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
phasei
label phasei
Definition: pEqn.H:27
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
PhaseCompressibleTurbulenceModel.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::fixedMultiPhaseHeatFluxFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fixedMultiPhaseHeatFluxFvPatchScalarField.C:114
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
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::fixedMultiPhaseHeatFluxFvPatchScalarField
Calculates a wall temperature that produces the specified overall wall heat flux across all the phase...
Definition: fixedMultiPhaseHeatFluxFvPatchScalarField.H:59
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::phase::name
const word & name() const
Definition: phase.H:111
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
fixedMultiPhaseHeatFluxFvPatchScalarField.H
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
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::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::fixedMultiPhaseHeatFluxFvPatchScalarField::fixedMultiPhaseHeatFluxFvPatchScalarField
fixedMultiPhaseHeatFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: fixedMultiPhaseHeatFluxFvPatchScalarField.C:41
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::fixedMultiPhaseHeatFluxFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: fixedMultiPhaseHeatFluxFvPatchScalarField.C:172
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