ENH: Added new heat transfer coefficient function object

Computes the heat transfer coefficient [W/m2/K] using a run-time
selectable model:
- ReynoldsAnalogy
- fixedReferenceTemperature
- localReferenceTemperature
This commit is contained in:
Andrew Heather
2017-12-22 19:00:25 +00:00
parent b85a38dc41
commit 38d13f41f3
11 changed files with 1611 additions and 0 deletions

View File

@ -0,0 +1,114 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "heatTransferCoeff.H"
#include "dictionary.H"
#include "heatTransferCoeffModel.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace functionObjects
{
defineTypeNameAndDebug(heatTransferCoeff, 0);
addToRunTimeSelectionTable(functionObject, heatTransferCoeff, dictionary);
}
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
bool Foam::functionObjects::heatTransferCoeff::calc()
{
volScalarField& htc = mesh_.lookupObjectRef<volScalarField>(resultName_);
htcModelPtr_->calc(htc);
return true;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::functionObjects::heatTransferCoeff::heatTransferCoeff
(
const word& name,
const Time& runTime,
const dictionary& dict
)
:
fieldExpression(name, runTime, dict),
htcModelPtr_()
{
read(dict);
setResultName(typeName, name + ":htc:" + htcModelPtr_->type());
volScalarField* heatTransferCoeffPtr
(
new volScalarField
(
IOobject
(
resultName_,
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh_,
dimensionedScalar("0", dimPower/dimArea/dimTemperature, 0.0)
)
);
mesh_.objectRegistry::store(heatTransferCoeffPtr);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::functionObjects::heatTransferCoeff::~heatTransferCoeff()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::functionObjects::heatTransferCoeff::read(const dictionary& dict)
{
if (fieldExpression::read(dict))
{
htcModelPtr_ = heatTransferCoeffModel::New(dict, mesh_, fieldName_);
htcModelPtr_->read(dict);
return true;
}
return false;
}
// ************************************************************************* //

View File

@ -0,0 +1,199 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
Class
Foam::functionObjects::heatTransferCoeff
Group
grpFieldFunctionObjects
Description
This function object calculates and writes the heat transfer coefficient
as a volScalarField for a set of patches.
The field is stored on the mesh database so that it can be retrieved and
used for other applications. Heat transfer coefficient, htc [W/m2/K]
can be evaluated using one of the following modes:
- ReynoldsAnalogy: Reynold's analogy
- localReferenceTemperature: local reference temperature
- fixedReferenceTemperature: specified reference temperature
Usage
Example usage for mode 'ReynoldsAnalogy' for incompressible case
\verbatim
htc
{
type heatTransferCoeff;
libs ("libfieldFunctionObjects.so");
field T;
patches ("walls.*");
htcModel ReynoldsAnalogy;
UInf (20 0 0);
Cp CpInf;
CpInf 1000;
rho rhoInf;
rhoInf 1.2;
}
\endverbatim
Example usage for mode 'ReynoldsAnalogy' for compressible case
\verbatim
htc
{
type heatTransferCoeff;
libs ("libfieldFunctionObjects.so");
field T;
patches ("walls.*");
htcModel ReynoldsAnalogy;
UInf (20 0 0);
}
\endverbatim
Example usage for mode 'localReferenceTemperature' for compressible case
\verbatim
htc
{
type heatTransferCoeff;
libs ("libfieldFunctionObjects.so");
field T;
patches ("walls.*");
htcModel local;
}
\endverbatim
Example usage for mode 'fixedReferenceTemperature' for compressible case
\verbatim
htc
{
type heatTransferCoeff;
libs ("libfieldFunctionObjects.so");
field T;
patches ("walls.*");
htcModel local;
TRef 300;
}
\endverbatim
See also
Foam::functionObjects::fieldExpression
Foam::heatTransferCoeffModels::fixedReferenceTemperature
Foam::heatTransferCoeffModels::localReferenceTemperature
SourceFiles
heatTransferCoeff.C
\*---------------------------------------------------------------------------*/
#ifndef functionObjects_heatTransferCoeff_H
#define functionObjects_heatTransferCoeff_H
#include "fieldExpression.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class heatTransferCoeffModel;
namespace functionObjects
{
/*---------------------------------------------------------------------------*\
Class heatTransferCoeff Declaration
\*---------------------------------------------------------------------------*/
class heatTransferCoeff
:
public fieldExpression
{
private:
// Private data
//- Heat transfer coefficient model
autoPtr<heatTransferCoeffModel> htcModelPtr_;
// Private Member Functions
//- Disallow default bitwise copy construct
heatTransferCoeff(const heatTransferCoeff&) = delete;
//- Disallow default bitwise assignment
void operator=(const heatTransferCoeff&) = delete;
protected:
//- Calculate the heat transfer coefficient field and return true
// if successful
virtual bool calc();
public:
//- Runtime type information
TypeName("heatTransferCoeff");
// Constructors
//- Construct for given objectRegistry and dictionary.
// Allow the possibility to load fields from files
heatTransferCoeff
(
const word& name,
const Time& runTime,
const dictionary& dict
);
//- Destructor
virtual ~heatTransferCoeff();
// Member Functions
//- Read the heatTransferCoeff data
virtual bool read(const dictionary&);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace functionObjects
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,271 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "ReynoldsAnalogy.H"
#include "fluidThermo.H"
#include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace heatTransferCoeffModels
{
defineTypeNameAndDebug(ReynoldsAnalogy, 0);
addToRunTimeSelectionTable
(
heatTransferCoeffModel,
ReynoldsAnalogy,
dictionary
);
}
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::Field<Foam::scalar>>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::rho(const label patchi) const
{
if (rhoName_ == "rhoInf")
{
const label n = mesh_.boundary()[patchi].size();
return tmp<Field<scalar>>(new Field<scalar>(n, rhoRef_));
}
else if (mesh_.foundObject<volScalarField>(rhoName_, false))
{
const volScalarField& rho =
mesh_.lookupObject<volScalarField>(rhoName_);
return rho.boundaryField()[patchi];
}
else
{
FatalErrorInFunction
<< "Unable to set rho for patch " << patchi
<< exit(FatalError);
}
return tmp<Field<scalar>>(nullptr);
}
Foam::tmp<Foam::Field<Foam::scalar>>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::Cp(const label patchi) const
{
if (CpName_ == "CpInf")
{
const label n = mesh_.boundary()[patchi].size();
return tmp<Field<scalar>>(new Field<scalar>(n, CpRef_));
}
else if (mesh_.foundObject<fluidThermo>(fluidThermo::typeName))
{
const fluidThermo& thermo =
mesh_.lookupObject<fluidThermo>(fluidThermo::typeName);
const scalarField& pp = thermo.p().boundaryField()[patchi];
const scalarField& Tp = thermo.T().boundaryField()[patchi];
return thermo.Cp(pp, Tp, patchi);
}
else
{
FatalErrorInFunction
<< "Unable to set Cp for patch " << patchi
<< exit(FatalError);
}
return tmp<Field<scalar>>(nullptr);
}
Foam::tmp<Foam::volSymmTensorField>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::devReff() const
{
typedef compressible::turbulenceModel cmpTurbModel;
typedef incompressible::turbulenceModel icoTurbModel;
if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
{
const cmpTurbModel& turb =
mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
return turb.devRhoReff()/turb.rho();
}
else if (mesh_.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
{
const incompressible::turbulenceModel& turb =
mesh_.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
return turb.devReff();
}
else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
{
const fluidThermo& thermo =
mesh_.lookupObject<fluidThermo>(fluidThermo::dictName);
const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
return -thermo.nu()*dev(twoSymm(fvc::grad(U)));
}
else if (mesh_.foundObject<transportModel>("transportProperties"))
{
const transportModel& laminarT =
mesh_.lookupObject<transportModel>("transportProperties");
const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
return -laminarT.nu()*dev(twoSymm(fvc::grad(U)));
}
else if (mesh_.foundObject<dictionary>("transportProperties"))
{
const dictionary& transportProperties =
mesh_.lookupObject<dictionary>("transportProperties");
dimensionedScalar nu(transportProperties.lookup("nu"));
const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
return -nu*dev(twoSymm(fvc::grad(U)));
}
else
{
FatalErrorInFunction
<< "No valid model for viscous stress calculation"
<< exit(FatalError);
return volSymmTensorField::null();
}
}
Foam::tmp<Foam::FieldField<Foam::Field, Foam::scalar>>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::Cf() const
{
const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
const volVectorField::Boundary& Ubf = U.boundaryField();
tmp<FieldField<Field, scalar>> tCf
(
new FieldField<Field, scalar>(Ubf.size())
);
FieldField<Field, scalar>& Cf = tCf.ref();
forAll(Cf, patchi)
{
Cf.set(patchi, new Field<scalar>(Ubf[patchi].size(), 0));
}
const volSymmTensorField R(devReff());
const volSymmTensorField::Boundary& Rbf = R.boundaryField();
for (label patchi : patchSet_)
{
const fvPatchVectorField& Up = Ubf[patchi];
const symmTensorField& Rp = Rbf[patchi];
const vectorField nHat(Up.patch().nf());
const scalarField tauByRhop(mag(nHat & Rp));
Cf[patchi] = 2*tauByRhop/magSqr(URef_);
}
return tCf;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::heatTransferCoeffModels::ReynoldsAnalogy::ReynoldsAnalogy
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
)
:
heatTransferCoeffModel(dict, mesh, TName),
UName_("U"),
URef_(vector::zero),
rhoName_("rho"),
rhoRef_(0.0),
CpName_("Cp"),
CpRef_(0.0)
{
read(dict);
}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
bool Foam::heatTransferCoeffModels::ReynoldsAnalogy::read
(
const dictionary& dict
)
{
if (heatTransferCoeffModel::read(dict))
{
dict.lookup("UInf") >> URef_;
dict.readIfPresent("Cp", CpName_);
if (CpName_ == "CpInf")
{
dict.lookup("CpInf") >> CpRef_;
}
dict.readIfPresent("rho", rhoName_);
if (rhoName_ == "rhoInf")
{
dict.lookup("rhoInf") >> rhoRef_;
}
return true;
}
return false;
}
void Foam::heatTransferCoeffModels::ReynoldsAnalogy::htc(volScalarField& htc)
{
const FieldField<Field, scalar> CfBf(Cf());
const scalar magU = mag(URef_);
volScalarField::Boundary& htcBf = htc.boundaryFieldRef();
forAllConstIters(patchSet_, iter)
{
label patchi = iter.key();
const scalarField rhop(rho(patchi));
const scalarField Cpp(Cp(patchi));
htcBf[patchi] = 0.5*rhop*Cpp*magU*CfBf[patchi];
}
}
// ************************************************************************* //

View File

@ -0,0 +1,182 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
Class
Foam::heatTransferCoeffModels::ReynoldsAnalogy
Description
Heat transfer coefficient calculation based on Reynolds Analogy
The heat transfer coefficient is derived from the skin friction
coefficient:
\f[
C_f = \frac{\tau_w}{0.5 \rho_\infty |U|^2}
\f]
as:
\f[
h = 0.5 \rho_\infty \C_{p,\infty} |U_{\infty}| C_f
\f]
Usage
Example of function object specification:
\verbatim
type heatTransferCoeff;
libs ("libfieldFunctionObjects.so");
...
htcModel ReynoldsAnalogy;
UInf (20 0 0);
Cp CpInf;
CpInf 1005;
...
\endverbatim
Where the entries comprise:
\table
Property | Description | Required | Default value
type | type name: heatTransferCoeff | yes |
htcModel | selected htc model | yes |
UInf | reference velocity | yes |
Cp | specific heat capacity field name | no |
rho | density field name | no |
\endtable
Note:
- to use a reference \c Cp, set \c Cp to \c CpInf
- to use a reference \c rho, set \c rho to \c rhoInf
SourceFiles
ReynoldsAnalogy.C
SeeAlso
Foam::heatTransferCoeffModel
\*---------------------------------------------------------------------------*/
#ifndef heatTransferCoeffModels_ReynoldsAnalogy_H
#define heatTransferCoeffModels_ReynoldsAnalogy_H
#include "heatTransferCoeffModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace heatTransferCoeffModels
{
/*---------------------------------------------------------------------------*\
Class ReynoldsAnalogy Declaration
\*---------------------------------------------------------------------------*/
class ReynoldsAnalogy
:
public heatTransferCoeffModel
{
// Private Member Functions
//- Disallow copy construct
ReynoldsAnalogy(const ReynoldsAnalogy&) = delete;
//- Disallow default bitwise assignment
void operator=(const ReynoldsAnalogy&) = delete;
protected:
// Protected data
//- Name of velocity field
word UName_;
//- Reference velocity
vector URef_;
//- Name of density field
word rhoName_;
//- Reference density
scalar rhoRef_;
//- Name of specific heat capacity field
word CpName_;
//- Reference specific heat capacity
scalar CpRef_;
// Protected Member Functions
virtual tmp<Field<scalar>> rho(const label patchi) const;
virtual tmp<Field<scalar>> Cp(const label patchi) const;
virtual tmp<volSymmTensorField> devReff() const;
tmp<FieldField<Field, scalar>> Cf() const;
//- Set the heat transfer coefficient
virtual void htc(volScalarField& htc);
public:
//- Runtime type information
TypeName("ReynoldsAnalogy");
// Constructors
//- Construct from components
ReynoldsAnalogy
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
);
//- Destructor
virtual ~ReynoldsAnalogy()
{}
// Member Functions
//- Read from dictionary
virtual bool read(const dictionary& dict);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace heatTransferCoeffModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,99 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "fixedReferenceTemperature.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace heatTransferCoeffModels
{
defineTypeNameAndDebug(fixedReferenceTemperature, 0);
addToRunTimeSelectionTable
(
heatTransferCoeffModel,
fixedReferenceTemperature,
dictionary
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::heatTransferCoeffModels::fixedReferenceTemperature::fixedReferenceTemperature
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
)
:
heatTransferCoeffModel(dict, mesh, TName),
TRef_(0)
{
read(dict);
}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
bool Foam::heatTransferCoeffModels::fixedReferenceTemperature::read
(
const dictionary& dict
)
{
if (heatTransferCoeffModel::read(dict))
{
dict.lookup("TRef") >> TRef_;
return true;
}
return false;
}
void Foam::heatTransferCoeffModels::fixedReferenceTemperature::htc
(
volScalarField& htc
)
{
const FieldField<Field, scalar> qBf(q());
const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
const volScalarField::Boundary& Tbf = T.boundaryField();
const scalar eps = ROOTVSMALL;
volScalarField::Boundary& htcBf = htc.boundaryFieldRef();
forAllConstIters(patchSet_, iter)
{
label patchi = iter.key();
htcBf[patchi] = qBf[patchi]/(TRef_ - Tbf[patchi] + eps);
}
}
// ************************************************************************* //

View File

@ -0,0 +1,145 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
Class
Foam::heatTransferCoeffModels::fixedReferenceTemperature
Description
Heat transfer coefficient calculation that employs a fixed reference
temperature
The heat transfer coefficient is specified by:
\f[
h = \frac{q}{T_{ref} - T_w}
\f]
Usage
Example of function object specification:
\verbatim
type heatTransferCoeff;
libs ("libfieldFunctionObjects.so");
...
htcModel fixedReferenceTemperature;
TRef 300;
...
\endverbatim
Where the entries comprise:
\table
Property | Description | Required | Default value
type | type name: heatTransferCoeff | yes |
htcModel | selected htc model | yes |
TRef | reference temperature | yes |
\endtable
SourceFiles
fixedReferenceTemperature.C
SeeAlso
Foam::heatTransferCoeffModel
\*---------------------------------------------------------------------------*/
#ifndef heatTransferCoeffModels_fixedReferenceTemperature_H
#define heatTransferCoeffModels_fixedReferenceTemperature_H
#include "heatTransferCoeffModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace heatTransferCoeffModels
{
/*---------------------------------------------------------------------------*\
Class fixedReferenceTemperature Declaration
\*---------------------------------------------------------------------------*/
class fixedReferenceTemperature
:
public heatTransferCoeffModel
{
// Private Member Functions
//- Disallow copy construct
fixedReferenceTemperature(const fixedReferenceTemperature&) = delete;
//- Disallow default bitwise assignment
void operator=(const fixedReferenceTemperature&) = delete;
protected:
// Protected data
//- Reference tempetaure
scalar TRef_;
// Protected Member Functions
//- Set the heat transfer coefficient
virtual void htc(volScalarField& htc);
public:
//- Runtime type information
TypeName("fixedReferenceTemperature");
// Constructors
//- Construct from components
fixedReferenceTemperature
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
);
//- Destructor
virtual ~fixedReferenceTemperature()
{}
// Member Functions
//- Read from dictionary
virtual bool read(const dictionary& dict);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace heatTransferCoeffModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,154 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "heatTransferCoeffModel.H"
#include "fvMesh.H"
#include "fluidThermo.H"
#include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(heatTransferCoeffModel, 0);
defineRunTimeSelectionTable(heatTransferCoeffModel, dictionary);
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::FieldField<Foam::Field, Foam::scalar>>
Foam::heatTransferCoeffModel::q() const
{
const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
const volScalarField::Boundary& Tbf = T.boundaryField();
tmp<FieldField<Field, scalar>> tq
(
new FieldField<Field, scalar>(Tbf.size())
);
FieldField<Field, scalar>& q = tq.ref();
forAll(q, patchi)
{
q.set(patchi, new Field<scalar>(Tbf[patchi].size(), 0));
}
typedef compressible::turbulenceModel cmpTurbModel;
if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
{
const cmpTurbModel& turb =
mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
const volScalarField& he = turb.transport().he();
const volScalarField::Boundary& hebf = he.boundaryField();
const volScalarField alphaEff(turb.alphaEff());
const volScalarField::Boundary& alphaEffbf = alphaEff.boundaryField();
for (label patchi : patchSet_)
{
q[patchi] = alphaEffbf[patchi]*hebf[patchi].snGrad();
}
}
else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
{
const fluidThermo& thermo =
mesh_.lookupObject<fluidThermo>(fluidThermo::dictName);
const volScalarField& he = thermo.he();
const volScalarField::Boundary& hebf = he.boundaryField();
const volScalarField& alpha(thermo.alpha());
const volScalarField::Boundary& alphabf = alpha.boundaryField();
for (label patchi : patchSet_)
{
q[patchi] = alphabf[patchi]*hebf[patchi].snGrad();
}
}
else
{
FatalErrorInFunction
<< "Unable to find a valid thermo model to evaluate q"
<< exit(FatalError);
}
// Add radiative heat flux contribution if present
if (mesh_.foundObject<volScalarField>(qrName_))
{
const volScalarField& qr = mesh_.lookupObject<volScalarField>(qrName_);
const volScalarField::Boundary& qrbf = qr.boundaryField();
for (label patchi : patchSet_)
{
q[patchi] += qrbf[patchi];
}
}
return tq;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::heatTransferCoeffModel::heatTransferCoeffModel
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
)
:
mesh_(mesh),
TName_(TName),
patchSet_(),
qrName_("qr")
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
bool Foam::heatTransferCoeffModel::read(const dictionary& dict)
{
const wordReList patchNames(dict.lookup("patches"));
patchSet_ = mesh_.boundaryMesh().patchSet(patchNames);
dict.readIfPresent("qr", qrName_);
return true;
}
bool Foam::heatTransferCoeffModel::calc(volScalarField& result)
{
htc(result);
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,158 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
Namespace
Foam::heatTransferCoeffModels
Description
A namespace for various heat transfer coefficient model implementations.
Class
Foam::heatTransferCoeffModel
Description
An abstract base class for heat transfer coeffcient models.
SourceFiles
heatTransferCoeffModel.C
heatTransferCoeffModelNew.C
\*---------------------------------------------------------------------------*/
#ifndef heatTransferCoeffModel_H
#define heatTransferCoeffModel_H
#include "dictionary.H"
#include "HashSet.H"
#include "volFields.H"
#include "runTimeSelectionTables.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class fvMesh;
/*---------------------------------------------------------------------------*\
Class heatTransferCoeffModel Declaration
\*---------------------------------------------------------------------------*/
class heatTransferCoeffModel
{
// Private Member Functions
//- Disallow copy construct
heatTransferCoeffModel(const heatTransferCoeffModel&) = delete;
//- Disallow default bitwise assignment
void operator=(const heatTransferCoeffModel&) = delete;
protected:
// Protected data
const fvMesh& mesh_;
const word TName_;
labelHashSet patchSet_;
word qrName_;
tmp<FieldField<Field, scalar>> q() const;
//- Set the heat transfer coefficient
virtual void htc(volScalarField& htc) = 0;
public:
//- Runtime type information
TypeName("heatTransferCoeffModel");
// Declare run-time constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
heatTransferCoeffModel,
dictionary,
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
),
(dict, mesh, TName)
);
// Selectors
//- Return a reference to the selected heat transfer coefficicent model
static autoPtr<heatTransferCoeffModel> New
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
);
// Constructors
//- Construct from components
heatTransferCoeffModel
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
);
//- Destructor
virtual ~heatTransferCoeffModel()
{}
// Member Functions
//- Read from dictionary
virtual bool read(const dictionary& dict);
virtual bool calc(volScalarField& result);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,59 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "heatTransferCoeffModel.H"
#include "fvMesh.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Foam::autoPtr<Foam::heatTransferCoeffModel> Foam::heatTransferCoeffModel::New
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
)
{
const word modelType(dict.lookup("htcModel"));
Info<< "Selecting heat transfer coefficient model " << modelType << endl;
auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
if (!cstrIter.found())
{
FatalErrorInFunction
<< "Unknown heatTransferCoeffModel type "
<< modelType << nl << nl
<< "Valid heatTransferCoeffModels :" << endl
<< dictionaryConstructorTablePtr_->sortedToc()
<< exit(FatalError);
}
return autoPtr<heatTransferCoeffModel>(cstrIter()(dict, mesh, TName));
}
// ************************************************************************* //

View File

@ -0,0 +1,93 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "localReferenceTemperature.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace heatTransferCoeffModels
{
defineTypeNameAndDebug(localReferenceTemperature, 0);
addToRunTimeSelectionTable
(
heatTransferCoeffModel,
localReferenceTemperature,
dictionary
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::heatTransferCoeffModels::localReferenceTemperature::
localReferenceTemperature
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
)
:
heatTransferCoeffModel(dict, mesh, TName)
{
read(dict);
}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
bool Foam::heatTransferCoeffModels::localReferenceTemperature::read
(
const dictionary& dict
)
{
return heatTransferCoeffModel::read(dict);
}
void Foam::heatTransferCoeffModels::localReferenceTemperature::htc
(
volScalarField& htc
)
{
const FieldField<Field, scalar> qBf(q());
const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
const volScalarField::Boundary& Tbf = T.boundaryField();
const scalar eps = ROOTVSMALL;
volScalarField::Boundary& htcBf = htc.boundaryFieldRef();
forAllConstIters(patchSet_, iter)
{
label patchi = iter.key();
const scalarField Tc(Tbf[patchi].patchInternalField());
htcBf[patchi] = qBf[patchi]/(Tc - Tbf[patchi] + eps);
}
}
// ************************************************************************* //

View File

@ -0,0 +1,137 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ 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 <http://www.gnu.org/licenses/>.
Class
Foam::heatTransferCoeffModels::localReferenceTemperature
Description
Heat transfer coefficient calculation that employs the patch internal
field as the reference temperature
The heat transfer coefficient is specified by:
\f[
h = \frac{q}{T_c - T_w}
\f]
Usage
Example of function object specification:
\verbatim
type heatTransferCoeff;
libs ("libfieldFunctionObjects.so");
...
htcModel localReferenceTemperature;
...
\endverbatim
Where the entries comprise:
\table
Property | Description | Required | Default value
type | type name: heatTransferCoeff | yes |
htcModel | selected htc model | yes |
\endtable
SourceFiles
localReferenceTemperature.C
SeeAlso
Foam::heatTransferCoeffModel
\*---------------------------------------------------------------------------*/
#ifndef heatTransferCoeffModels_localReferenceTemperature_H
#define heatTransferCoeffModels_localReferenceTemperature_H
#include "heatTransferCoeffModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace heatTransferCoeffModels
{
/*---------------------------------------------------------------------------*\
Class localReferenceTemperature Declaration
\*---------------------------------------------------------------------------*/
class localReferenceTemperature
:
public heatTransferCoeffModel
{
// Private Member Functions
//- Disallow copy construct
localReferenceTemperature(const localReferenceTemperature&) = delete;
//- Disallow default bitwise assignment
void operator=(const localReferenceTemperature&) = delete;
protected:
// Protected Member Functions
//- Set the heat transfer coefficient
virtual void htc(volScalarField& htc);
public:
//- Runtime type information
TypeName("localReferenceTemperature");
// Constructors
//- Construct from components
localReferenceTemperature
(
const dictionary& dict,
const fvMesh& mesh,
const word& TName
);
//- Destructor
virtual ~localReferenceTemperature()
{}
// Member Functions
//- Read from dictionary
virtual bool read(const dictionary& dict);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace heatTransferCoeffModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //