/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see .
\*---------------------------------------------------------------------------*/
#include "physicoChemicalConstants.H"
using namespace Foam::constant;
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template
inline const Foam::ThermoCloud&
Foam::ThermoCloud::cloudCopy() const
{
return cloudCopyPtr_();
}
template
inline const typename CloudType::particleType::constantProperties&
Foam::ThermoCloud::constProps() const
{
return constProps_;
}
template
inline typename CloudType::particleType::constantProperties&
Foam::ThermoCloud::constProps()
{
return constProps_;
}
template
inline const Foam::SLGThermo& Foam::ThermoCloud::thermo() const
{
return thermo_;
}
template
inline const Foam::volScalarField& Foam::ThermoCloud::T() const
{
return T_;
}
template
inline const Foam::volScalarField& Foam::ThermoCloud::p() const
{
return p_;
}
template
inline const Foam::HeatTransferModel >&
Foam::ThermoCloud::heatTransfer() const
{
return heatTransferModel_;
}
template
inline const Foam::scalarIntegrationScheme&
Foam::ThermoCloud::TIntegrator() const
{
return TIntegrator_;
}
template
inline bool Foam::ThermoCloud::radiation() const
{
return radiation_;
}
template
inline Foam::DimensionedField&
Foam::ThermoCloud::radAreaP()
{
if (!radiation_)
{
FatalErrorIn
(
"inline Foam::DimensionedField "
"Foam::ThermoCloud::radAreaP()"
) << "Radiation field requested, but radiation model not active"
<< abort(FatalError);
}
return radAreaP_();
}
template
inline const Foam::DimensionedField&
Foam::ThermoCloud::radAreaP() const
{
if (!radiation_)
{
FatalErrorIn
(
"inline Foam::DimensionedField "
"Foam::ThermoCloud::radAreaP()"
) << "Radiation field requested, but radiation model not active"
<< abort(FatalError);
}
return radAreaP_();
}
template
inline Foam::DimensionedField&
Foam::ThermoCloud::radT4()
{
if (!radiation_)
{
FatalErrorIn
(
"inline Foam::DimensionedField "
"Foam::ThermoCloud::radT4()"
) << "Radiation field requested, but radiation model not active"
<< abort(FatalError);
}
return radT4_();
}
template
inline const Foam::DimensionedField&
Foam::ThermoCloud::radT4() const
{
if (!radiation_)
{
FatalErrorIn
(
"inline Foam::DimensionedField "
"Foam::ThermoCloud::radT4()"
) << "Radiation field requested, but radiation model not active"
<< abort(FatalError);
}
return radT4_();
}
template
inline Foam::DimensionedField&
Foam::ThermoCloud::radAreaPT4()
{
if (!radiation_)
{
FatalErrorIn
(
"inline Foam::DimensionedField "
"Foam::ThermoCloud::radAreaPT4()"
) << "Radiation field requested, but radiation model not active"
<< abort(FatalError);
}
return radAreaPT4_();
}
template
inline const Foam::DimensionedField&
Foam::ThermoCloud::radAreaPT4() const
{
if (!radiation_)
{
FatalErrorIn
(
"inline Foam::DimensionedField "
"Foam::ThermoCloud::radAreaPT4()"
) << "Radiation field requested, but radiation model not active"
<< abort(FatalError);
}
return radAreaPT4_();
}
template
inline Foam::DimensionedField&
Foam::ThermoCloud::hsTrans()
{
return hsTrans_();
}
template
inline const Foam::DimensionedField&
Foam::ThermoCloud::hsTrans() const
{
return hsTrans_();
}
template
inline Foam::DimensionedField&
Foam::ThermoCloud::hsCoeff()
{
return hsCoeff_();
}
template
inline const Foam::DimensionedField&
Foam::ThermoCloud::hsCoeff() const
{
return hsCoeff_();
}
template
inline Foam::tmp
Foam::ThermoCloud::Sh(volScalarField& hs) const
{
if (debug)
{
Info<< "hsTrans min/max = " << min(hsTrans()).value() << ", "
<< max(hsTrans()).value() << nl
<< "hsCoeff min/max = " << min(hsCoeff()).value() << ", "
<< max(hsCoeff()).value() << endl;
}
if (this->solution().coupled())
{
if (this->solution().semiImplicit("h"))
{
const volScalarField Cp(thermo_.thermo().Cp());
const DimensionedField
Vdt(this->mesh().V()*this->db().time().deltaT());
return
hsTrans()/Vdt
- fvm::SuSp(hsCoeff()/(Cp*Vdt), hs)
+ hsCoeff()/(Cp*Vdt)*hs;
}
else
{
tmp tfvm(new fvScalarMatrix(hs, dimEnergy/dimTime));
fvScalarMatrix& fvm = tfvm();
fvm.source() = -hsTrans()/(this->db().time().deltaT());
return tfvm;
}
}
return tmp(new fvScalarMatrix(hs, dimEnergy/dimTime));
}
template
inline Foam::tmp Foam::ThermoCloud::Ep() const
{
tmp tEp
(
new volScalarField
(
IOobject
(
this->name() + ":radiation:Ep",
this->db().time().timeName(),
this->db(),
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh(),
dimensionedScalar("zero", dimMass/dimLength/pow3(dimTime), 0.0)
)
);
if (radiation_)
{
scalarField& Ep = tEp().internalField();
const scalar dt = this->db().time().deltaTValue();
const scalarField& V = this->mesh().V();
const scalar epsilon = constProps_.epsilon0();
const scalarField& sumAreaPT4 = radAreaPT4_->field();
Ep = sumAreaPT4*epsilon*physicoChemical::sigma.value()/V/dt;
}
return tEp;
}
template
inline Foam::tmp Foam::ThermoCloud::ap() const
{
tmp tap
(
new volScalarField
(
IOobject
(
this->name() + ":radiation:ap",
this->db().time().timeName(),
this->db(),
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh(),
dimensionedScalar("zero", dimless/dimLength, 0.0)
)
);
if (radiation_)
{
scalarField& ap = tap().internalField();
const scalar dt = this->db().time().deltaTValue();
const scalarField& V = this->mesh().V();
const scalar epsilon = constProps_.epsilon0();
const scalarField& sumAreaP = radAreaP_->field();
ap = sumAreaP*epsilon/V/dt;
}
return tap;
}
template
inline Foam::tmp
Foam::ThermoCloud::sigmap() const
{
tmp tsigmap
(
new volScalarField
(
IOobject
(
this->name() + ":radiation:sigmap",
this->db().time().timeName(),
this->db(),
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh(),
dimensionedScalar("zero", dimless/dimLength, 0.0)
)
);
if (radiation_)
{
scalarField& sigmap = tsigmap().internalField();
const scalar dt = this->db().time().deltaTValue();
const scalarField& V = this->mesh().V();
const scalar epsilon = constProps_.epsilon0();
const scalar f = constProps_.f0();
const scalarField& sumAreaP = radAreaP_->field();
sigmap *= sumAreaP*(1.0 - f)*(1.0 - epsilon)/V/dt;
}
return tsigmap;
}
template
inline Foam::scalar Foam::ThermoCloud::Tmax() const
{
scalar T = -GREAT;
scalar n = 0;
forAllConstIter(typename ThermoCloud, *this, iter)
{
const parcelType& p = iter();
T = max(T, p.T());
n++;
}
reduce(T, maxOp());
reduce(n, sumOp