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-2018 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"
32 #include "ThermalDiffusivity.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
41  const DimensionedField<scalar, volMesh>& iF
42 )
43 :
44  fixedValueFvPatchScalarField(p, iF),
45  q_(p.size(), 0.0),
46  relax_(1.0),
47  Tmin_(0.0)
48 {}
49 
50 
53 (
54  const fvPatch& p,
55  const DimensionedField<scalar, volMesh>& iF,
56  const dictionary& dict
57 )
58 :
59  fixedValueFvPatchScalarField(p, iF, dict),
60  q_("q", dict, p.size()),
61  relax_(dict.lookupOrDefault<scalar>("relax", 1.0)),
62  Tmin_(dict.lookupOrDefault<scalar>("Tmin", 273))
63 {}
64 
65 
68 (
69  const fixedMultiPhaseHeatFluxFvPatchScalarField& psf,
70  const fvPatch& p,
71  const DimensionedField<scalar, volMesh>& iF,
72  const fvPatchFieldMapper& mapper
73 )
74 :
75  fixedValueFvPatchScalarField(psf, p, iF, mapper),
76  q_(psf.q_, mapper),
77  relax_(psf.relax_),
78  Tmin_(psf.Tmin_)
79 {}
80 
81 
84 (
85  const fixedMultiPhaseHeatFluxFvPatchScalarField& psf
86 )
87 :
88  fixedValueFvPatchScalarField(psf),
89  q_(psf.q_),
90  relax_(psf.relax_),
91  Tmin_(psf.Tmin_)
92 {}
93 
94 
97 (
98  const fixedMultiPhaseHeatFluxFvPatchScalarField& psf,
99  const DimensionedField<scalar, volMesh>& iF
100 )
101 :
102  fixedValueFvPatchScalarField(psf, iF),
103  q_(psf.q_),
104  relax_(psf.relax_),
105  Tmin_(psf.Tmin_)
106 {}
107 
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 {
114  if (updated())
115  {
116  return;
117  }
118 
119  // Lookup the fluid model
120  const phaseSystem& fluid =
121  (
122  db().lookupObject<phaseSystem>("phaseProperties")
123  );
124 
125  const scalarField& Tp = *this;
126 
127  scalarField A(Tp.size(), scalar(0));
128  scalarField B(Tp.size(), scalar(0));
129  scalarField Q(Tp.size(), scalar(0));
130 
131  forAll(fluid.phases(), phasei)
132  {
133  const phaseModel& phase = fluid.phases()[phasei];
134  const fluidThermo& thermo = phase.thermo();
135 
136  const fvPatchScalarField& alpha =
137  phase.boundaryField()[patch().index()];
138 
139  const fvPatchScalarField& T =
140  thermo.T().boundaryField()[patch().index()];
141 
142  const scalarField kappaEff(phase.kappaEff(patch().index()));
143 
144  if (debug)
145  {
146  scalarField q0(T.snGrad()*alpha*kappaEff);
147  Q += q0;
148 
149  Info<< patch().name() << " " << phase.name()
150  << ": Heat flux " << gMin(q0) << " - " << gMax(q0) << endl;
151  }
152 
153  A += T.patchInternalField()*alpha*kappaEff*patch().deltaCoeffs();
154  B += alpha*kappaEff*patch().deltaCoeffs();
155  }
156 
157  if (debug)
158  {
159  Info<< patch().name() << " " << ": overall heat flux "
160  << gMin(Q) << " - " << gMax(Q) << " W/m2, power: "
161  << gSum(patch().magSf()*Q) << " W" << endl;
162  }
163 
164  operator==((1 - relax_)*Tp + relax_*max(Tmin_,(q_ + A)/(B)));
165 
166  fixedValueFvPatchScalarField::updateCoeffs();
167 }
168 
169 
171 {
173  os.writeKeyword("relax") << relax_ << token::END_STATEMENT << nl;
174  q_.writeEntry("q", os);
175  writeEntry("value", os);
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
186  fixedMultiPhaseHeatFluxFvPatchScalarField
187  );
188 }
189 
190 
191 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
label phasei
Definition: pEqn.H:27
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMin(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static const char nl
Definition: Ostream.H:265
Type gMax(const FieldField< Field, Type > &f)
messageStream Info
volScalarField & p
virtual void write(Ostream &) const
Write.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.