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-2020 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 )
86 :
87  fixedValueFvPatchScalarField(psf),
88  q_(psf.q_),
89  relax_(psf.relax_),
90  Tmin_(psf.Tmin_)
91 {}
92 
93 
96 (
97  const fixedMultiPhaseHeatFluxFvPatchScalarField& psf,
98  const DimensionedField<scalar, volMesh>& iF
99 )
100 :
101  fixedValueFvPatchScalarField(psf, iF),
102  q_(psf.q_),
103  relax_(psf.relax_),
104  Tmin_(psf.Tmin_)
105 {}
106 
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 {
113  if (updated())
114  {
115  return;
116  }
117 
118  // Lookup the fluid model
119  const phaseSystem& fluid =
120  (
121  db().lookupObject<phaseSystem>("phaseProperties")
122  );
123 
124  const scalarField& Tp = *this;
125 
126  scalarField A(Tp.size(), scalar(0));
127  scalarField B(Tp.size(), scalar(0));
128  scalarField Q(Tp.size(), scalar(0));
129 
130  forAll(fluid.phases(), phasei)
131  {
132  const phaseModel& phase = fluid.phases()[phasei];
133  const fluidThermo& thermo = phase.thermo();
134 
135  const fvPatchScalarField& alpha =
136  phase.boundaryField()[patch().index()];
137 
138  const fvPatchScalarField& T =
139  thermo.T().boundaryField()[patch().index()];
140 
141  const scalarField kappaEff(phase.kappaEff(patch().index()));
142 
143  if (debug)
144  {
145  scalarField q0(T.snGrad()*alpha*kappaEff);
146  Q += q0;
147 
148  Info<< patch().name() << " " << phase.name()
149  << ": Heat flux " << gMin(q0) << " - " << gMax(q0) << endl;
150  }
151 
152  A += T.patchInternalField()*alpha*kappaEff*patch().deltaCoeffs();
153  B += alpha*kappaEff*patch().deltaCoeffs();
154  }
155 
156  if (debug)
157  {
158  Info<< patch().name() << " " << ": overall heat flux "
159  << gMin(Q) << " - " << gMax(Q) << " W/m2, power: "
160  << gSum(patch().magSf()*Q) << " W" << endl;
161  }
162 
163  operator==((1 - relax_)*Tp + relax_*max(Tmin_,(q_ + A)/(B)));
164 
165  fixedValueFvPatchScalarField::updateCoeffs();
166 }
167 
168 
170 (
171  const fvPatchFieldMapper& m
172 )
173 {
174  fixedValueFvPatchScalarField::autoMap(m);
175  m(q_, q_);
176 }
177 
178 
180 (
181  const fvPatchScalarField& ptf,
182  const labelList& addr
183 )
184 {
185  fixedValueFvPatchScalarField::rmap(ptf, addr);
186 
187  const fixedMultiPhaseHeatFluxFvPatchScalarField& mptf =
188  refCast<const fixedMultiPhaseHeatFluxFvPatchScalarField>(ptf);
189 
190  q_.rmap(mptf.q_, addr);
191 }
192 
193 
195 {
197  writeEntry(os, "relax", relax_);
198  writeEntry(os, "q", q_);
199  writeEntry(os, "value", *this);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 namespace Foam
206 {
208  (
210  fixedMultiPhaseHeatFluxFvPatchScalarField
211  );
212 }
213 
214 
215 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
label phasei
Definition: pEqn.H:27
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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:280
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.
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.