fixedMultiPhaseHeatFluxFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2022 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 
30 #include "phaseSystem.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
40  const DimensionedField<scalar, volMesh>& iF
41 )
42 :
43  fixedValueFvPatchScalarField(p, iF),
44  q_(p.size(), 0.0),
45  relax_(1.0),
46  Tmin_(0.0)
47 {}
48 
49 
52 (
53  const fvPatch& p,
54  const DimensionedField<scalar, volMesh>& iF,
55  const dictionary& dict
56 )
57 :
58  fixedValueFvPatchScalarField(p, iF, dict),
59  q_("q", dict, p.size()),
60  relax_(dict.lookupOrDefault<scalar>("relax", 1.0)),
61  Tmin_(dict.lookupOrDefault<scalar>("Tmin", 273))
62 {}
63 
64 
67 (
68  const fixedMultiPhaseHeatFluxFvPatchScalarField& psf,
69  const fvPatch& p,
70  const DimensionedField<scalar, volMesh>& iF,
71  const fvPatchFieldMapper& mapper
72 )
73 :
74  fixedValueFvPatchScalarField(psf, p, iF, mapper),
75  q_(mapper(psf.q_)),
76  relax_(psf.relax_),
77  Tmin_(psf.Tmin_)
78 {}
79 
80 
83 (
84  const fixedMultiPhaseHeatFluxFvPatchScalarField& psf,
85  const DimensionedField<scalar, volMesh>& iF
86 )
87 :
88  fixedValueFvPatchScalarField(psf, iF),
89  q_(psf.q_),
90  relax_(psf.relax_),
91  Tmin_(psf.Tmin_)
92 {}
93 
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
100  if (updated())
101  {
102  return;
103  }
104 
105  // Lookup the fluid model
106  const phaseSystem& fluid =
107  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
108 
109  const scalarField& Tp = *this;
110 
111  scalarField A(Tp.size(), scalar(0));
112  scalarField B(Tp.size(), scalar(0));
113  scalarField Q(Tp.size(), scalar(0));
114 
115  forAll(fluid.phases(), phasei)
116  {
117  const phaseModel& phase = fluid.phases()[phasei];
118  const fluidThermo& thermo = phase.thermo();
119 
120  const fvPatchScalarField& alpha =
121  phase.boundaryField()[patch().index()];
122 
123  const fvPatchScalarField& T =
124  thermo.T().boundaryField()[patch().index()];
125 
126  const scalarField kappaEff(phase.kappaEff(patch().index()));
127 
128  if (debug)
129  {
130  scalarField q0(T.snGrad()*alpha*kappaEff);
131  Q += q0;
132 
133  Info<< patch().name() << " " << phase.name()
134  << ": Heat flux " << gMin(q0) << " - " << gMax(q0) << endl;
135  }
136 
137  A += T.patchInternalField()*alpha*kappaEff*patch().deltaCoeffs();
138  B += alpha*kappaEff*patch().deltaCoeffs();
139  }
140 
141  if (debug)
142  {
143  Info<< patch().name() << " " << ": overall heat flux "
144  << gMin(Q) << " - " << gMax(Q) << " W/m2, power: "
145  << gSum(patch().magSf()*Q) << " W" << endl;
146  }
147 
148  operator==((1 - relax_)*Tp + relax_*max(Tmin_,(q_ + A)/(B)));
149 
150  fixedValueFvPatchScalarField::updateCoeffs();
151 }
152 
153 
155 (
156  const fvPatchFieldMapper& m
157 )
158 {
159  fixedValueFvPatchScalarField::autoMap(m);
160  m(q_, q_);
161 }
162 
163 
165 (
166  const fvPatchScalarField& ptf,
167  const labelList& addr
168 )
169 {
170  fixedValueFvPatchScalarField::rmap(ptf, addr);
171 
172  const fixedMultiPhaseHeatFluxFvPatchScalarField& mptf =
173  refCast<const fixedMultiPhaseHeatFluxFvPatchScalarField>(ptf);
174 
175  q_.rmap(mptf.q_, addr);
176 }
177 
178 
180 (
181  const fvPatchScalarField& ptf
182 )
183 {
184  fixedValueFvPatchScalarField::reset(ptf);
185 
186  const fixedMultiPhaseHeatFluxFvPatchScalarField& mptf =
187  refCast<const fixedMultiPhaseHeatFluxFvPatchScalarField>(ptf);
188 
189  q_.reset(mptf.q_);
190 }
191 
192 
194 {
196  writeEntry(os, "relax", relax_);
197  writeEntry(os, "q", q_);
198  writeEntry(os, "value", *this);
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 namespace Foam
205 {
207  (
209  fixedMultiPhaseHeatFluxFvPatchScalarField
210  );
211 }
212 
213 
214 // ************************************************************************* //
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:237
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
label phasei
Definition: pEqn.H:27
Type gMin(const FieldField< Field, Type > &f)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
Type gSum(const FieldField< Field, Type > &f)
fixedMultiPhaseHeatFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
fvPatchField< scalar > fvPatchScalarField
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
messageStream Info
volScalarField & p
virtual void write(Ostream &) const
Write.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.