lumpedMassWallTemperatureFvPatchScalarField.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) 2016-2019 OpenCFD Ltd.
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 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "mappedPatchBase.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  mixedFvPatchScalarField(p, iF),
44  temperatureCoupledBase(patch()), // default method (fluidThermo)
45  Cp_(0.0),
46  mass_(0.0),
47  curTimeIndex_(-1)
48 {
49  refValue() = 0.0;
50  refGrad() = 0.0;
51  valueFraction() = 1.0;
52 }
53 
54 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  mixedFvPatchScalarField(ptf, p, iF, mapper),
66  Cp_(ptf.Cp_),
67  mass_(ptf.mass_),
68  curTimeIndex_(-1)
69 {}
70 
71 
74 (
75  const fvPatch& p,
77  const dictionary& dict
78 )
79 :
80  mixedFvPatchScalarField(p, iF),
82  Cp_(dict.get<scalar>("Cp")),
83  mass_(dict.get<scalar>("mass")),
84  curTimeIndex_(-1)
85 {
87  this->readValueEntry(dict, IOobjectOption::MUST_READ);
88  refValue() = *this;
89  refGrad() = Zero;
90  valueFraction() = 1.0;
91 }
92 
93 
96 (
98 )
99 :
100  mixedFvPatchScalarField(tppsf),
101  temperatureCoupledBase(tppsf),
102  Cp_(tppsf.Cp_),
103  mass_(tppsf.mass_),
104  curTimeIndex_(-1)
105 {}
106 
107 
110 (
113 )
114 :
115  mixedFvPatchScalarField(tppsf, iF),
116  temperatureCoupledBase(patch(), tppsf),
117  Cp_(tppsf.Cp_),
118  mass_(tppsf.mass_),
119  curTimeIndex_(-1)
120 {}
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 (
127  const fvPatchFieldMapper& mapper
128 )
129 {
130  mixedFvPatchScalarField::autoMap(mapper);
132 }
133 
134 
136 (
137  const fvPatchField<scalar>& ptf,
138  const labelList& addr
139 )
140 {
141  mixedFvPatchScalarField::rmap(ptf, addr);
142 
144  refCast
145  <
147  >(ptf);
148 
149  temperatureCoupledBase::rmap(tiptf, addr);
150 }
151 
152 
154 {
155  if (updated() || (curTimeIndex_ == this->db().time().timeIndex()))
156  {
157  return;
158  }
159 
160  // Calculate heat flux in or out the wall
161  scalarField& Tp(*this);
162  const scalarField& magSf = patch().magSf();
163 
164  const scalar deltaT(db().time().deltaTValue());
165 
166  tmp<scalarField> tkappa(kappa(Tp));
167 
168  const scalarField q(tkappa.ref()*snGrad());
169 
170  // Total heat in or out of the wall
171  const scalar Q = gSum(q*magSf);
172 
173  Tp += -(Q/mass_/Cp_)*deltaT;
174 
175  refGrad() = 0.0;
176  refValue() = Tp;
177  valueFraction() = 1.0;
178 
179  mixedFvPatchScalarField::updateCoeffs();
180 
181  if (debug)
182  {
183  scalar Qin(0);
184  scalar Qout(0);
185 
186  forAll(q, facei)
187  {
188  if (q[facei] > 0.0) // out the wall
189  {
190  Qout += q[facei]*magSf[facei];
191  }
192  else if (q[facei] < 0.0) // into the wall
193  {
194  Qin += q[facei]*magSf[facei];
195  }
196  }
197 
198  Info<< patch().boundaryMesh().mesh().name() << ':'
199  << patch().name() << ':'
200  << this->internalField().name() << " :"
201  << " heat transfer rate:" << Q
202  << " wall temperature "
203  << " min:" << gMin(*this)
204  << " max:" << gMax(*this)
205  << " avg:" << gAverage(*this)
206  << " Qin [W]:" << Qin
207  << " Qout [W]:" << Qout
208  << endl;
209  }
210 
211  curTimeIndex_ = this->db().time().timeIndex();
212 }
213 
214 
216 (
217  Ostream& os
218 ) const
219 {
222 
223  os.writeEntry("Cp", Cp_);
224  os.writeEntry("mass", mass_);
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 namespace Foam
231 {
233  (
235  lumpedMassWallTemperatureFvPatchScalarField
236  );
237 }
238 
239 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Type gMin(const FieldField< Field, Type > &f)
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
virtual void readDict(const dictionary &dict)
Read dictionary entries.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Type gSum(const FieldField< Field, Type > &f)
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
Common functions used in temperature coupled boundaries.
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual void write(Ostream &) const
Write.
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
lumpedMassWallTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.