Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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)
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),
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
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 = gWeightedSum(magSf, q);
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 auto limits = gMinMax(*this);
199 auto avg = gAverage(*this);
200
201 Info<< patch().boundaryMesh().mesh().name() << ':'
202 << patch().name() << ':'
203 << this->internalField().name() << " :"
204 << " heat transfer rate:" << Q
205 << " wall temperature "
206 << " min:" << limits.min()
207 << " max:" << limits.max()
208 << " avg:" << avg
209 << " Qin [W]:" << Qin
210 << " Qout [W]:" << Qout
212 }
213
214 curTimeIndex_ = this->db().time().timeIndex();
215}
216
217
219(
220 Ostream& os
221) const
222{
225
226 os.writeEntry("Cp", Cp_);
227 os.writeEntry("mass", mass_);
229
230
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233namespace Foam
234{
236 (
238 lumpedMassWallTemperatureFvPatchScalarField
239 );
240}
241
242// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void readDict(const dictionary &dict)
Read dictionary entries.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
lumpedMassWallTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
Common functions used in temperature coupled boundaries.
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 rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
temperatureCoupledBase(const fvPatch &patch, const KMethodType method=KMethodType::mtFluidThermo)
Default construct from patch, using fluidThermo (default) or specified method.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
volScalarField & p
auto limits
Definition setRDeltaT.H:186
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< scalar > fvPatchScalarField
label timeIndex
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299