ThermoParcel.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) 2011-2017 OpenFOAM Foundation
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 
28 #include "ThermoParcel.H"
30 
31 using namespace Foam::constant;
32 
33 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
34 
35 template<class ParcelType>
36 template<class TrackCloudType>
38 (
39  TrackCloudType& cloud,
40  trackingData& td
41 )
42 {
43  ParcelType::setCellValues(cloud, td);
44 
45  tetIndices tetIs = this->currentTetIndices();
46 
47  td.Cpc() = td.CpInterp().interpolate(this->coordinates(), tetIs);
48 
49  td.Tc() = td.TInterp().interpolate(this->coordinates(), tetIs);
50 
51  if (td.Tc() < cloud.constProps().TMin())
52  {
53  if (debug)
54  {
56  << "Limiting observed temperature in cell " << this->cell()
57  << " to " << cloud.constProps().TMin() << nl << endl;
58  }
59 
60  td.Tc() = cloud.constProps().TMin();
61  }
62 }
63 
64 
65 template<class ParcelType>
66 template<class TrackCloudType>
68 (
69  TrackCloudType& cloud,
70  trackingData& td,
71  const scalar dt
72 )
73 {
74  const label celli = this->cell();
75  const scalar massCell = this->massCell(td);
76 
77  td.Uc() += cloud.UTrans()[celli]/massCell;
78 
79 // tetIndices tetIs = this->currentTetIndices();
80 // Tc_ = td.TInterp().interpolate(this->coordinates(), tetIs);
81 
82  const scalar CpMean = td.CpInterp().psi()[celli];
83  td.Tc() += cloud.hsTrans()[celli]/(CpMean*massCell);
84 
85  if (td.Tc() < cloud.constProps().TMin())
86  {
87  if (debug)
88  {
90  << "Limiting observed temperature in cell " << celli
91  << " to " << cloud.constProps().TMin() << nl << endl;
92  }
93 
94  td.Tc() = cloud.constProps().TMin();
95  }
96 }
97 
98 
99 template<class ParcelType>
100 template<class TrackCloudType>
102 (
103  TrackCloudType& cloud,
104  trackingData& td,
105  const scalar T,
106  scalar& Ts,
107  scalar& rhos,
108  scalar& mus,
109  scalar& Pr,
110  scalar& kappas
111 ) const
112 {
113  // Surface temperature using two thirds rule
114  Ts = (2.0*T + td.Tc())/3.0;
115 
116  if (Ts < cloud.constProps().TMin())
117  {
118  if (debug)
119  {
121  << "Limiting parcel surface temperature to "
122  << cloud.constProps().TMin() << nl << endl;
123  }
124 
125  Ts = cloud.constProps().TMin();
126  }
127 
128  // Assuming thermo props vary linearly with T for small d(T)
129  const scalar TRatio = td.Tc()/Ts;
130 
131  rhos = td.rhoc()*TRatio;
132 
133  tetIndices tetIs = this->currentTetIndices();
134  mus = td.muInterp().interpolate(this->coordinates(), tetIs)/TRatio;
135  kappas = td.kappaInterp().interpolate(this->coordinates(), tetIs)/TRatio;
136 
137  Pr = td.Cpc()*mus/kappas;
138  Pr = max(ROOTVSMALL, Pr);
139 }
140 
141 
142 template<class ParcelType>
143 template<class TrackCloudType>
145 (
146  TrackCloudType& cloud,
147  trackingData& td,
148  const scalar dt
149 )
150 {
151  // Define local properties at beginning of time step
152  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
153  const scalar np0 = this->nParticle_;
154  const scalar mass0 = this->mass();
155 
156  // Store T for consistent radiation source
157  const scalar T0 = this->T_;
158 
159 
160  // Calc surface values
161  // ~~~~~~~~~~~~~~~~~~~
162  scalar Ts, rhos, mus, Pr, kappas;
163  calcSurfaceValues(cloud, td, this->T_, Ts, rhos, mus, Pr, kappas);
164 
165  // Reynolds number
166  scalar Re = this->Re(rhos, this->U_, td.Uc(), this->d_, mus);
167 
168 
169  // Sources
170  // ~~~~~~~
171 
172  // Explicit momentum source for particle
173  vector Su = Zero;
174 
175  // Linearised momentum source coefficient
176  scalar Spu = 0.0;
177 
178  // Momentum transfer from the particle to the carrier phase
179  vector dUTrans = Zero;
180 
181  // Explicit enthalpy source for particle
182  scalar Sh = 0.0;
183 
184  // Linearised enthalpy source coefficient
185  scalar Sph = 0.0;
186 
187  // Sensible enthalpy transfer from the particle to the carrier phase
188  scalar dhsTrans = 0.0;
189 
190 
191  // Heat transfer
192  // ~~~~~~~~~~~~~
193 
194  // Sum Ni*Cpi*Wi of emission species
195  scalar NCpW = 0.0;
196 
197  // Calculate new particle temperature
198  this->T_ =
199  this->calcHeatTransfer
200  (
201  cloud,
202  td,
203  dt,
204  Re,
205  Pr,
206  kappas,
207  NCpW,
208  Sh,
209  dhsTrans,
210  Sph
211  );
212 
213 
214  // Motion
215  // ~~~~~~
216 
217  // Calculate new particle velocity
218  this->U_ =
219  this->calcVelocity(cloud, td, dt, Re, mus, mass0, Su, dUTrans, Spu);
220 
221 
222  // Accumulate carrier phase source terms
223  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
224  if (cloud.solution().coupled())
225  {
226  // Update momentum transfer
227  cloud.UTrans()[this->cell()] += np0*dUTrans;
228 
229  // Update momentum transfer coefficient
230  cloud.UCoeff()[this->cell()] += np0*Spu;
231 
232  // Update sensible enthalpy transfer
233  cloud.hsTrans()[this->cell()] += np0*dhsTrans;
234 
235  // Update sensible enthalpy coefficient
236  cloud.hsCoeff()[this->cell()] += np0*Sph;
237 
238  // Update radiation fields
239  if (cloud.radiation())
240  {
241  const scalar ap = this->areaP();
242  const scalar T4 = pow4(T0);
243  cloud.radAreaP()[this->cell()] += dt*np0*ap;
244  cloud.radT4()[this->cell()] += dt*np0*T4;
245  cloud.radAreaPT4()[this->cell()] += dt*np0*ap*T4;
246  }
247  }
248 }
249 
250 
251 template<class ParcelType>
252 template<class TrackCloudType>
254 (
255  TrackCloudType& cloud,
256  trackingData& td,
257  const scalar dt,
258  const scalar Re,
259  const scalar Pr,
260  const scalar kappa,
261  const scalar NCpW,
262  const scalar Sh,
263  scalar& dhsTrans,
264  scalar& Sph
265 )
266 {
267  if (!cloud.heatTransfer().active())
268  {
269  return T_;
270  }
271 
272  const scalar d = this->d();
273  const scalar rho = this->rho();
274  const scalar As = this->areaS(d);
275  const scalar V = this->volume(d);
276  const scalar m = rho*V;
277 
278  // Calc heat transfer coefficient
279  scalar htc = cloud.heatTransfer().htc(d, Re, Pr, kappa, NCpW);
280 
281  // Calculate the integration coefficients
282  const scalar bcp = htc*As/(m*Cp_);
283  const scalar acp = bcp*td.Tc();
284  scalar ancp = Sh;
285  if (cloud.radiation())
286  {
287  const tetIndices tetIs = this->currentTetIndices();
288  const scalar Gc = td.GInterp().interpolate(this->coordinates(), tetIs);
289  const scalar sigma = physicoChemical::sigma.value();
290  const scalar epsilon = cloud.constProps().epsilon0();
291 
292  ancp += As*epsilon*(Gc/4.0 - sigma*pow4(T_));
293  }
294  ancp /= m*Cp_;
295 
296  // Integrate to find the new parcel temperature
297  const scalar deltaT = cloud.TIntegrator().delta(T_, dt, acp + ancp, bcp);
298  const scalar deltaTncp = ancp*dt;
299  const scalar deltaTcp = deltaT - deltaTncp;
300 
301  // Calculate the new temperature and the enthalpy transfer terms
302  scalar Tnew = T_ + deltaT;
303  Tnew = min(max(Tnew, cloud.constProps().TMin()), cloud.constProps().TMax());
304 
305  dhsTrans -= m*Cp_*deltaTcp;
306 
307  Sph = dt*m*Cp_*bcp;
308 
309  return Tnew;
310 }
311 
312 
313 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
314 
315 template<class ParcelType>
317 (
319 )
320 :
321  ParcelType(p),
322  T_(p.T_),
323  Cp_(p.Cp_)
324 {}
325 
326 
327 template<class ParcelType>
329 (
330  const ThermoParcel<ParcelType>& p,
331  const polyMesh& mesh
332 )
333 :
334  ParcelType(p, mesh),
335  T_(p.T_),
336  Cp_(p.Cp_)
337 {}
338 
339 
340 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
341 
342 #include "ThermoParcelIO.C"
343 
344 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::ThermoParcel::trackingData::kappaInterp
const interpolation< scalar > & kappaInterp() const
Return const access to the interpolator for continuous.
Definition: ThermoParcelTrackingDataI.H:117
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::ThermoParcel::trackingData::TInterp
const interpolation< scalar > & TInterp() const
Return const access to the interpolator for continuous.
Definition: ThermoParcelTrackingDataI.H:101
Foam::ThermoParcel::trackingData::CpInterp
const interpolation< scalar > & CpInterp() const
Return const access to the interpolator for continuous.
Definition: ThermoParcelTrackingDataI.H:109
Foam::ThermoParcel::trackingData::GInterp
const interpolation< scalar > & GInterp() const
Return const access to the interpolator for continuous.
Definition: ThermoParcelTrackingDataI.H:125
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::ThermoParcel::calcHeatTransfer
scalar calcHeatTransfer(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW, const scalar Sh, scalar &dhsTrans, scalar &Sph)
Calculate new particle temperature.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fieldTypes::volume
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
rho
rho
Definition: readInitialConditions.H:88
Foam::ThermoParcel::calc
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Definition: ThermoParcel.C:145
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::ThermoParcel::calcSurfaceValues
void calcSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
Definition: ThermoParcel.C:102
Foam::interpolation::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Return the field to be interpolated.
Definition: interpolation.H:127
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::ThermoParcel::trackingData
Definition: ThermoParcel.H:152
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::ThermoParcel::setCellValues
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Definition: ThermoParcel.C:38
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
ThermoParcelIO.C
Foam::ThermoParcel::cellValueSourceCorrection
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Definition: ThermoParcel.C:68
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::ThermoParcel
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
Definition: ThermoParcel.H:57
physicoChemicalConstants.H
Foam::Vector< scalar >
Foam::ThermoParcel::trackingData::Cpc
scalar Cpc() const
Return the continuous phase specific heat capacity.
Definition: ThermoParcelTrackingDataI.H:153
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::ThermoParcel::trackingData::Tc
scalar Tc() const
Return the continuous phase temperature.
Definition: ThermoParcelTrackingDataI.H:139
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
ThermoParcel.H
Foam::interpolation::interpolate
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
T0
scalar T0
Definition: createFields.H:22
Foam::ThermoParcel::ThermoParcel
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: ThermoParcelI.H:77
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54