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 | Copyright (C) 2011-2015 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 
26 #include "ThermoParcel.H"
28 
29 using namespace Foam::constant;
30 
31 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
32 
33 template<class ParcelType>
34 template<class TrackData>
36 (
37  TrackData& td,
38  const scalar dt,
39  const label cellI
40 )
41 {
42  ParcelType::setCellValues(td, dt, cellI);
43 
44  tetIndices tetIs = this->currentTetIndices();
45 
46  Cpc_ = td.CpInterp().interpolate(this->position(), tetIs);
47 
48  Tc_ = td.TInterp().interpolate(this->position(), tetIs);
49 
50  if (Tc_ < td.cloud().constProps().TMin())
51  {
52  if (debug)
53  {
55  << "Limiting observed temperature in cell " << cellI << " to "
56  << td.cloud().constProps().TMin() << nl << endl;
57  }
58 
59  Tc_ = td.cloud().constProps().TMin();
60  }
61 }
62 
63 
64 template<class ParcelType>
65 template<class TrackData>
67 (
68  TrackData& td,
69  const scalar dt,
70  const label cellI
71 )
72 {
73  this->Uc_ += td.cloud().UTrans()[cellI]/this->massCell(cellI);
74 
75  const scalar CpMean = td.CpInterp().psi()[cellI];
76 
77  tetIndices tetIs = this->currentTetIndices();
78  Tc_ = td.TInterp().interpolate(this->position(), tetIs);
79  Tc_ += td.cloud().hsTrans()[cellI]/(CpMean*this->massCell(cellI));
80 
81  if (Tc_ < td.cloud().constProps().TMin())
82  {
83  if (debug)
84  {
86  << "Limiting observed temperature in cell " << cellI << " to "
87  << td.cloud().constProps().TMin() << nl << endl;
88  }
89 
90  Tc_ = td.cloud().constProps().TMin();
91  }
92 }
93 
94 
95 template<class ParcelType>
96 template<class TrackData>
98 (
99  TrackData& td,
100  const label cellI,
101  const scalar T,
102  scalar& Ts,
103  scalar& rhos,
104  scalar& mus,
105  scalar& Pr,
106  scalar& kappas
107 ) const
108 {
109  // Surface temperature using two thirds rule
110  Ts = (2.0*T + Tc_)/3.0;
111 
112  if (Ts < td.cloud().constProps().TMin())
113  {
114  if (debug)
115  {
117  << "Limiting parcel surface temperature to "
118  << td.cloud().constProps().TMin() << nl << endl;
119  }
120 
121  Ts = td.cloud().constProps().TMin();
122  }
123 
124  // Assuming thermo props vary linearly with T for small d(T)
125  const scalar TRatio = Tc_/Ts;
126 
127  rhos = this->rhoc_*TRatio;
128 
129  tetIndices tetIs = this->currentTetIndices();
130  mus = td.muInterp().interpolate(this->position(), tetIs)/TRatio;
131  kappas = td.kappaInterp().interpolate(this->position(), tetIs)/TRatio;
132 
133  Pr = Cpc_*mus/kappas;
134  Pr = max(ROOTVSMALL, Pr);
135 }
136 
137 
138 template<class ParcelType>
139 template<class TrackData>
141 (
142  TrackData& td,
143  const scalar dt,
144  const label cellI
145 )
146 {
147  // Define local properties at beginning of time step
148  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
149  const scalar np0 = this->nParticle_;
150  const scalar mass0 = this->mass();
151 
152  // Store T for consistent radiation source
153  const scalar T0 = this->T_;
154 
155 
156  // Calc surface values
157  // ~~~~~~~~~~~~~~~~~~~
158  scalar Ts, rhos, mus, Pr, kappas;
159  calcSurfaceValues(td, cellI, this->T_, Ts, rhos, mus, Pr, kappas);
160 
161  // Reynolds number
162  scalar Re = this->Re(this->U_, this->d_, rhos, mus);
163 
164 
165  // Sources
166  // ~~~~~~~
167 
168  // Explicit momentum source for particle
170 
171  // Linearised momentum source coefficient
172  scalar Spu = 0.0;
173 
174  // Momentum transfer from the particle to the carrier phase
175  vector dUTrans = vector::zero;
176 
177  // Explicit enthalpy source for particle
178  scalar Sh = 0.0;
179 
180  // Linearised enthalpy source coefficient
181  scalar Sph = 0.0;
182 
183  // Sensible enthalpy transfer from the particle to the carrier phase
184  scalar dhsTrans = 0.0;
185 
186 
187  // Heat transfer
188  // ~~~~~~~~~~~~~
189 
190  // Sum Ni*Cpi*Wi of emission species
191  scalar NCpW = 0.0;
192 
193  // Calculate new particle temperature
194  this->T_ =
195  this->calcHeatTransfer
196  (
197  td,
198  dt,
199  cellI,
200  Re,
201  Pr,
202  kappas,
203  NCpW,
204  Sh,
205  dhsTrans,
206  Sph
207  );
208 
209 
210  // Motion
211  // ~~~~~~
212 
213  // Calculate new particle velocity
214  this->U_ =
215  this->calcVelocity(td, dt, cellI, Re, mus, mass0, Su, dUTrans, Spu);
216 
217 
218  // Accumulate carrier phase source terms
219  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
220  if (td.cloud().solution().coupled())
221  {
222  // Update momentum transfer
223  td.cloud().UTrans()[cellI] += np0*dUTrans;
224 
225  // Update momentum transfer coefficient
226  td.cloud().UCoeff()[cellI] += np0*Spu;
227 
228  // Update sensible enthalpy transfer
229  td.cloud().hsTrans()[cellI] += np0*dhsTrans;
230 
231  // Update sensible enthalpy coefficient
232  td.cloud().hsCoeff()[cellI] += np0*Sph;
233 
234  // Update radiation fields
235  if (td.cloud().radiation())
236  {
237  const scalar ap = this->areaP();
238  const scalar T4 = pow4(T0);
239  td.cloud().radAreaP()[cellI] += dt*np0*ap;
240  td.cloud().radT4()[cellI] += dt*np0*T4;
241  td.cloud().radAreaPT4()[cellI] += dt*np0*ap*T4;
242  }
243  }
244 }
245 
246 
247 template<class ParcelType>
248 template<class TrackData>
250 (
251  TrackData& td,
252  const scalar dt,
253  const label cellI,
254  const scalar Re,
255  const scalar Pr,
256  const scalar kappa,
257  const scalar NCpW,
258  const scalar Sh,
259  scalar& dhsTrans,
260  scalar& Sph
261 )
262 {
263  if (!td.cloud().heatTransfer().active())
264  {
265  return T_;
266  }
267 
268  const scalar d = this->d();
269  const scalar rho = this->rho();
270 
271  // Calc heat transfer coefficient
272  scalar htc = td.cloud().heatTransfer().htc(d, Re, Pr, kappa, NCpW);
273 
274  if (mag(htc) < ROOTVSMALL && !td.cloud().radiation())
275  {
276  return
277  max
278  (
279  T_ + dt*Sh/(this->volume(d)*rho*Cp_),
280  td.cloud().constProps().TMin()
281  );
282  }
283 
284  htc = max(htc, ROOTVSMALL);
285  const scalar As = this->areaS(d);
286 
287  scalar ap = Tc_ + Sh/(As*htc);
288  scalar bp = 6.0*(Sh/As + htc*(Tc_ - T_));
289  if (td.cloud().radiation())
290  {
291  tetIndices tetIs = this->currentTetIndices();
292  const scalar Gc = td.GInterp().interpolate(this->position(), tetIs);
293  const scalar sigma = physicoChemical::sigma.value();
294  const scalar epsilon = td.cloud().constProps().epsilon0();
295 
296  // Assume constant source
297  scalar s = epsilon*(Gc/4.0 - sigma*pow4(T_));
298 
299  ap += s/htc;
300  bp += 6.0*s;
301  }
302  bp /= rho*d*Cp_*(ap - T_) + ROOTVSMALL;
303 
304  // Integrate to find the new parcel temperature
306  td.cloud().TIntegrator().integrate(T_, dt, ap*bp, bp);
307 
308  scalar Tnew =
309  min
310  (
311  max
312  (
313  Tres.value(),
314  td.cloud().constProps().TMin()
315  ),
316  td.cloud().constProps().TMax()
317  );
318 
319  Sph = dt*htc*As;
320 
321  dhsTrans += Sph*(Tres.average() - Tc_);
322 
323  return Tnew;
324 }
325 
326 
327 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
328 
329 template<class ParcelType>
331 (
333 )
334 :
335  ParcelType(p),
336  T_(p.T_),
337  Cp_(p.Cp_),
338  Tc_(p.Tc_),
339  Cpc_(p.Cpc_)
340 {}
341 
342 
343 template<class ParcelType>
345 (
346  const ThermoParcel<ParcelType>& p,
347  const polyMesh& mesh
348 )
349 :
350  ParcelType(p, mesh),
351  T_(p.T_),
352  Cp_(p.Cp_),
353  Tc_(p.Tc_),
354  Cpc_(p.Cpc_)
355 {}
356 
357 
358 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
359 
360 #include "ThermoParcelIO.C"
361 
362 // ************************************************************************* //
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
Sh
scalar Sh
Definition: solveChemistry.H:2
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
Foam::Re
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Foam::ThermoParcel::calc
void calc(TrackData &td, const scalar dt, const label cellI)
Update parcel properties over the time interval.
Definition: ThermoParcel.C:141
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::ThermoParcel::setCellValues
void setCellValues(TrackData &td, const scalar dt, const label cellI)
Set cell values.
Definition: ThermoParcel.C:36
Foam::ThermoParcel::calcSurfaceValues
void calcSurfaceValues(TrackData &td, const label cellI, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
Definition: ThermoParcel.C:98
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
ThermoParcelIO.C
Foam::IntegrationScheme::integrationResult
Helper class to supply results of integration.
Definition: IntegrationScheme.H:57
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::ThermoParcel::calcHeatTransfer
scalar calcHeatTransfer(TrackData &td, const scalar dt, const label cellI, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW, const scalar Sh, scalar &dhsTrans, scalar &Sph)
Calculate new particle temperature.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:73
epsilon
epsilon
Definition: createTDFields.H:56
Foam::ThermoParcel::ThermoParcel
ThermoParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
Definition: ThermoParcelI.H:75
Foam::ThermoParcel
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
Definition: ThermoParcel.H:51
physicoChemicalConstants.H
Foam::Vector< scalar >
Foam::ThermoParcel::cellValueSourceCorrection
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label cellI)
Correct cell values using latest transfer information.
Definition: ThermoParcel.C:67
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
ThermoParcel.H
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259