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),
45  (
46  patch(),
47  "undefined",
48  "undefined",
49  "undefined-K",
50  "undefined-alpha"
51  ),
52  Cp_(0.0),
53  mass_(0.0),
54  curTimeIndex_(-1)
55 {
56  refValue() = 0.0;
57  refGrad() = 0.0;
58  valueFraction() = 1.0;
59 }
60 
61 
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  mixedFvPatchScalarField(ptf, p, iF, mapper),
73  Cp_(ptf.Cp_),
74  mass_(ptf.mass_),
75  curTimeIndex_(-1)
76 {}
77 
78 
81 (
82  const fvPatch& p,
84  const dictionary& dict
85 )
86 :
87  mixedFvPatchScalarField(p, iF),
89  Cp_(dict.get<scalar>("Cp")),
90  mass_(dict.get<scalar>("mass")),
91  curTimeIndex_(-1)
92 {
93  refGrad() = 0.0;
94  valueFraction() = 1.0;
95  refValue() = scalarField("value", dict, p.size());
96 
97  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
98 }
99 
100 
103 (
105 )
106 :
107  mixedFvPatchScalarField(tppsf),
108  temperatureCoupledBase(tppsf),
109  Cp_(tppsf.Cp_),
110  mass_(tppsf.mass_),
111  curTimeIndex_(-1)
112 {}
113 
114 
117 (
120 )
121 :
122  mixedFvPatchScalarField(tppsf, iF),
123  temperatureCoupledBase(patch(), tppsf),
124  Cp_(tppsf.Cp_),
125  mass_(tppsf.mass_),
126  curTimeIndex_(-1)
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 (
134  const fvPatchFieldMapper& mapper
135 )
136 {
137  mixedFvPatchScalarField::autoMap(mapper);
138  temperatureCoupledBase::autoMap(mapper);
139 }
140 
141 
143 (
144  const fvPatchField<scalar>& ptf,
145  const labelList& addr
146 )
147 {
148  mixedFvPatchScalarField::rmap(ptf, addr);
149 
151  refCast
152  <
154  >(ptf);
155 
156  temperatureCoupledBase::rmap(tiptf, addr);
157 }
158 
159 
161 {
162  if (updated() || (curTimeIndex_ == this->db().time().timeIndex()))
163  {
164  return;
165  }
166 
167  // Calculate heat flux in or out the wall
168  scalarField& Tp(*this);
169  const scalarField& magSf = patch().magSf();
170 
171  const scalar deltaT(db().time().deltaTValue());
172 
173  tmp<scalarField> tkappa(kappa(Tp));
174 
175  const scalarField q(tkappa.ref()*snGrad());
176 
177  // Total heat in or out of the wall
178  const scalar Q = gSum(q*magSf);
179 
180  Tp += -(Q/mass_/Cp_)*deltaT;
181 
182  refGrad() = 0.0;
183  refValue() = Tp;
184  valueFraction() = 1.0;
185 
186  mixedFvPatchScalarField::updateCoeffs();
187 
188  if (debug)
189  {
190  scalar Qin(0);
191  scalar Qout(0);
192 
193  forAll(q, facei)
194  {
195  if (q[facei] > 0.0) // out the wall
196  {
197  Qout += q[facei]*magSf[facei];
198  }
199  else if (q[facei] < 0.0) // into the wall
200  {
201  Qin += q[facei]*magSf[facei];
202  }
203  }
204 
205  Info<< patch().boundaryMesh().mesh().name() << ':'
206  << patch().name() << ':'
207  << this->internalField().name() << " :"
208  << " heat transfer rate:" << Q
209  << " wall temperature "
210  << " min:" << gMin(*this)
211  << " max:" << gMax(*this)
212  << " avg:" << gAverage(*this)
213  << " Qin [W]:" << Qin
214  << " Qout [W]:" << Qout
215  << endl;
216  }
217 
218  curTimeIndex_ = this->db().time().timeIndex();
219 }
220 
221 
223 (
224  Ostream& os
225 ) const
226 {
229 
230  os.writeEntry("Cp", Cp_);
231  os.writeEntry("mass", mass_);
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 namespace Foam
238 {
240  (
243  );
244 }
245 
246 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
volFields.H
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
lumpedMassWallTemperatureFvPatchScalarField.H
Foam::lumpedMassWallTemperatureFvPatchScalarField::rmap
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: lumpedMassWallTemperatureFvPatchScalarField.C:143
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:135
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::lumpedMassWallTemperatureFvPatchScalarField
Employs a lumped mass model for temperature.
Definition: lumpedMassWallTemperatureFvPatchScalarField.H:98
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::lumpedMassWallTemperatureFvPatchScalarField::lumpedMassWallTemperatureFvPatchScalarField
lumpedMassWallTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: lumpedMassWallTemperatureFvPatchScalarField.C:38
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Foam::lumpedMassWallTemperatureFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: lumpedMassWallTemperatureFvPatchScalarField.C:133
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::temperatureCoupledBase::kappa
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Definition: temperatureCoupledBase.C:210
Foam::lumpedMassWallTemperatureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: lumpedMassWallTemperatureFvPatchScalarField.C:160
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
Foam::temperatureCoupledBase::write
void write(Ostream &os) const
Write.
Definition: temperatureCoupledBase.C:512
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::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::lumpedMassWallTemperatureFvPatchScalarField::write
void write(Ostream &) const
Write.
Definition: lumpedMassWallTemperatureFvPatchScalarField.C:223
mappedPatchBase.H