solidDisplacementFoam: Updated to solve for the solid temperature in energy conservative form

consistent with the solid solver module and solidThermo.
This commit is contained in:
Henry Weller
2022-10-26 16:31:25 +01:00
parent 714e13a970
commit f05bcb541c
9 changed files with 27 additions and 739 deletions

View File

@ -12,8 +12,6 @@ if (thermo.thermalStress())
const volScalarField& E(thermo.E()); const volScalarField& E(thermo.E());
const volScalarField& nu(thermo.nu()); const volScalarField& nu(thermo.nu());
const volScalarField Cp(thermo.Cp());
Info<< "Calculating Lame's coefficients\n" << endl; Info<< "Calculating Lame's coefficients\n" << endl;
const volScalarField mu(E/(2*(1 + nu))); const volScalarField mu(E/(2*(1 + nu)));

View File

@ -73,20 +73,25 @@ int main(int argc, char *argv[])
if (thermo.thermalStress()) if (thermo.thermalStress())
{ {
volScalarField& T = thermo.T(); volScalarField& e = thermo.he();
fvScalarMatrix TEqn
fvScalarMatrix eEqn
( (
fvm::ddt(rho, Cp, T) fvm::ddt(rho, e)
+ thermophysicalTransport->divq(T) + thermophysicalTransport->divq(e)
== ==
fvModels.source(rho*Cp, T) fvModels.source(rho, e)
); );
fvConstraints.constrain(TEqn); eEqn.relax();
TEqn.solve(); fvConstraints.constrain(eEqn);
fvConstraints.constrain(T); eEqn.solve();
fvConstraints.constrain(e);
thermo.correct();
} }
{ {

View File

@ -24,12 +24,8 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "solidDisplacementThermo.H" #include "solidDisplacementThermo.H"
#include "fvMesh.H"
#include "fvmLaplacian.H"
#include "fvcSnGrad.H"
#include "surfaceInterpolate.H"
/* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */ /* * * * * * * * * * * * * * * Private Static Data * * * * * * * * * * * * * */
namespace Foam namespace Foam
{ {
@ -37,45 +33,6 @@ namespace Foam
} }
void Foam::solidDisplacementThermo::readProperty(volScalarField& prop) const
{
const dictionary& propDict(subDict(prop.name()));
const word propType(propDict.lookup("type"));
if (propType == "uniform")
{
prop == dimensionedScalar
(
prop.name(),
prop.dimensions(),
propDict.lookup<scalar>("value")
);
}
else if (propType == "field")
{
const volScalarField propField
(
IOobject
(
prop.name(),
prop.mesh().time().timeName(0),
prop.mesh(),
IOobject::MUST_READ,
IOobject::NO_WRITE
),
prop.mesh()
);
prop == propField;
}
else
{
FatalErrorInFunction
<< "Valid type entries are uniform or field for " << prop.name()
<< abort(FatalError);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::solidDisplacementThermo::solidDisplacementThermo Foam::solidDisplacementThermo::solidDisplacementThermo
@ -84,82 +41,13 @@ Foam::solidDisplacementThermo::solidDisplacementThermo
const word& phaseName const word& phaseName
) )
: :
solidThermo::composite(mesh, phaseName), constSolidThermo(mesh, phaseName),
planeStress_(lookup("planeStress")), planeStress_(lookup("planeStress")),
thermalStress_(lookup("thermalStress")), thermalStress_(lookup("thermalStress")),
Cp_ E_(readProperty("E", dimPressure)),
( nu_(readProperty("nu", dimless)),
IOobject alphav_(readProperty("alphav", dimless/dimTemperature))
( {}
phasePropertyName("Cp"),
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimEnergy/dimMass/dimTemperature
),
kappa_
(
IOobject
(
phasePropertyName("kappa"),
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
Cp_.dimensions()*dimensionSet(1, -1, -1, 0, 0)
),
E_
(
IOobject
(
phasePropertyName("E"),
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimPressure
),
nu_
(
IOobject
(
phasePropertyName("nu"),
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimless
),
alphav_
(
IOobject
(
phasePropertyName("alphav"),
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimless/dimTemperature
)
{
readProperty(rho_);
readProperty(Cp_);
readProperty(kappa_);
readProperty(E_);
readProperty(nu_);
readProperty(alphav_);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
@ -170,21 +58,6 @@ Foam::solidDisplacementThermo::~solidDisplacementThermo()
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::rho() const
{
return rho_;
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::rho
(
const label patchi
) const
{
return rho_.boundaryField()[patchi];
}
const Foam::volScalarField& Foam::solidDisplacementThermo::E() const const Foam::volScalarField& Foam::solidDisplacementThermo::E() const
{ {
return E_; return E_;
@ -230,337 +103,4 @@ const Foam::scalarField& Foam::solidDisplacementThermo::alphav
} }
Foam::volScalarField& Foam::solidDisplacementThermo::he()
{
NotImplemented;
return rho_;
}
const Foam::volScalarField& Foam::solidDisplacementThermo::he() const
{
NotImplemented;
return rho_;
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::he
(
const scalarField& T,
const labelList& cells
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::he
(
const volScalarField& p,
const volScalarField& T
) const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::he
(
const scalarField& T,
const label patchi
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::hs() const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::hs
(
const volScalarField& p,
const volScalarField& T
) const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::hs
(
const scalarField& T,
const label patchi
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::hs
(
const scalarField& T,
const labelList& cells
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::ha() const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::ha
(
const volScalarField& p,
const volScalarField& T
) const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::ha
(
const scalarField& T,
const label patchi
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::ha
(
const scalarField& T,
const labelList& cells
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::hc() const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::THE
(
const volScalarField& h,
const volScalarField& p,
const volScalarField& T0
) const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::THE
(
const scalarField& he,
const scalarField& T0,
const labelList& cells
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::THE
(
const scalarField& he,
const scalarField& T0,
const label patchi
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
const Foam::volScalarField& Foam::solidDisplacementThermo::Cp() const
{
return Cp_;
}
const Foam::volScalarField& Foam::solidDisplacementThermo::Cv() const
{
return Cp_;
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::Cp
(
const scalarField& T,
const label patchi
) const
{
return Cp_.boundaryField()[patchi];
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::Cv
(
const scalarField& T,
const label patchi
) const
{
return Cp_.boundaryField()[patchi];
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::Cpv() const
{
return Cp_;
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::Cpv
(
const scalarField& T,
const label patchi
) const
{
return Cp_.boundaryField()[patchi];
}
const Foam::volScalarField& Foam::solidDisplacementThermo::kappa() const
{
return kappa_;
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::alphahe() const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::alphahe
(
const label patchi
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::kappaEff
(
const volScalarField& alphat
) const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::kappaEff
(
const scalarField& alphat,
const label patchi
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::volScalarField> Foam::solidDisplacementThermo::alphaEff
(
const volScalarField& alphat
) const
{
NotImplemented;
return tmp<volScalarField>(nullptr);
}
Foam::tmp<Foam::scalarField> Foam::solidDisplacementThermo::alphaEff
(
const scalarField& alphat,
const label patchi
) const
{
NotImplemented;
return tmp<scalarField>(nullptr);
}
Foam::tmp<Foam::volVectorField> Foam::solidDisplacementThermo::Kappa() const
{
NotImplemented;
return tmp<volVectorField>(nullptr);
}
Foam::tmp<Foam::vectorField> Foam::solidDisplacementThermo::Kappa
(
const label patchi
) const
{
NotImplemented;
return tmp<vectorField>(nullptr);
}
Foam::tmp<Foam::symmTensorField> Foam::solidDisplacementThermo::KappaLocal
(
const label patchi
) const
{
NotImplemented;
return tmp<symmTensorField>(nullptr);
}
void Foam::solidDisplacementThermo::correct()
{}
bool Foam::solidDisplacementThermo::read()
{
return regIOobject::read();
}
Foam::tmp<Foam::surfaceScalarField>
Foam::solidDisplacementThermo::q() const
{
return -fvc::interpolate(kappa_)*fvc::snGrad(T());
}
Foam::tmp<Foam::fvScalarMatrix>
Foam::solidDisplacementThermo::divq(volScalarField& T) const
{
return -fvm::laplacian(kappa_, T);
}
// ************************************************************************* // // ************************************************************************* //

View File

@ -35,7 +35,7 @@ SourceFiles
#ifndef solidDisplacementThermo_H #ifndef solidDisplacementThermo_H
#define solidDisplacementThermo_H #define solidDisplacementThermo_H
#include "solidThermo.H" #include "constSolidThermo.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -48,7 +48,7 @@ namespace Foam
class solidDisplacementThermo class solidDisplacementThermo
: :
public solidThermo::composite public constSolidThermo
{ {
// Private data // Private data
@ -58,12 +58,6 @@ class solidDisplacementThermo
//- Switch to enable thermal stress //- Switch to enable thermal stress
Switch thermalStress_; Switch thermalStress_;
//- Heat capacity at constant pressure [J/kg/K]
volScalarField Cp_;
//- Thermal conductivity [W/m/K]
volScalarField kappa_;
//- Youngs modulus [Pa] //- Youngs modulus [Pa]
volScalarField E_; volScalarField E_;
@ -74,11 +68,6 @@ class solidDisplacementThermo
volScalarField alphav_; volScalarField alphav_;
// Private Member Functions
void readProperty(volScalarField& prop) const;
public: public:
//- Runtime type information //- Runtime type information
@ -101,26 +90,6 @@ public:
// Member Functions // Member Functions
//- Return the name of the thermo physics
virtual word thermoName() const
{
return type();
}
//- Return true if the equation of state is incompressible
// i.e. rho != f(p)
virtual bool incompressible() const
{
return true;
}
//- Return true if the equation of state is isochoric
// i.e. rho = const
virtual bool isochoric() const
{
return true;
}
//- Returns true to enable plane stress //- Returns true to enable plane stress
bool planeStress() const bool planeStress() const
{ {
@ -136,12 +105,6 @@ public:
// Access to thermophysical state variables // Access to thermophysical state variables
//- Density [kg/m^3]
virtual tmp<volScalarField> rho() const;
//- Density for patch [kg/m^3]
virtual tmp<scalarField> rho(const label patchi) const;
//- Youngs modulus [Pa] //- Youngs modulus [Pa]
virtual const volScalarField& E() const; virtual const volScalarField& E() const;
@ -159,221 +122,6 @@ public:
//- Volumetric thermal expansion coefficient for a patch [1/T] //- Volumetric thermal expansion coefficient for a patch [1/T]
virtual const scalarField& alphav(const label patchi) const; virtual const scalarField& alphav(const label patchi) const;
//- Enthalpy/Internal energy [J/kg]
// Non-const access allowed for transport equations
virtual volScalarField& he();
//- Enthalpy/Internal energy [J/kg]
virtual const volScalarField& he() const;
//- Heat capacity at constant pressure [J/kg/K]
virtual const volScalarField& Cp() const;
//- Heat capacity at constant volume [J/kg/K]
virtual const volScalarField& Cv() const;
// Access to transport state variables
//- Thermal conductivity of mixture [W/m/K]
virtual const volScalarField& kappa() const;
// Fields derived from thermodynamic state variables
//- Enthalpy/Internal energy
// for given pressure and temperature [J/kg]
virtual tmp<volScalarField> he
(
const volScalarField& p,
const volScalarField& T
) const;
//- Enthalpy/Internal energy for cell-set [J/kg]
virtual tmp<scalarField> he
(
const scalarField& T,
const labelList& cells
) const;
//- Enthalpy/Internal energy for patch [J/kg]
virtual tmp<scalarField> he
(
const scalarField& T,
const label patchi
) const;
//- Sensible enthalpy [J/kg]
virtual tmp<volScalarField> hs() const;
//- Sensible enthalpy
// for given pressure and temperature [J/kg]
virtual tmp<volScalarField> hs
(
const volScalarField& p,
const volScalarField& T
) const;
//- Sensible enthalpy for cell-set [J/kg]
virtual tmp<scalarField> hs
(
const scalarField& T,
const labelList& cells
) const;
//- Sensible enthalpy for patch [J/kg]
virtual tmp<scalarField> hs
(
const scalarField& T,
const label patchi
) const;
//- Absolute enthalpy [J/kg/K]
virtual tmp<volScalarField> ha() const;
//- Absolute enthalpy
// for given pressure and temperature [J/kg]
virtual tmp<volScalarField> ha
(
const volScalarField& p,
const volScalarField& T
) const;
//- Absolute enthalpy for patch [J/kg/K]
virtual tmp<scalarField> ha
(
const scalarField& T,
const label patchi
) const;
//- Absolute enthalpy for cell-set [J/kg]
virtual tmp<scalarField> ha
(
const scalarField& T,
const labelList& cells
) const;
//- Enthalpy of formation [J/kg]
virtual tmp<volScalarField> hc() const;
//- Temperature from enthalpy/internal energy
virtual tmp<volScalarField> THE
(
const volScalarField& h,
const volScalarField& p,
const volScalarField& T0 // starting temperature
) const;
//- Temperature from enthalpy/internal energy for cell-set
virtual tmp<scalarField> THE
(
const scalarField& he,
const scalarField& T0, // starting temperature
const labelList& cells
) const;
//- Temperature from enthalpy/internal energy for patch
virtual tmp<scalarField> THE
(
const scalarField& he,
const scalarField& T0, // starting temperature
const label patchi
) const;
//- Heat capacity at constant pressure for patch [J/kg/K]
virtual tmp<scalarField> Cp
(
const scalarField& T,
const label patchi
) const;
//- Heat capacity at constant volume for patch [J/kg/K]
virtual tmp<scalarField> Cv
(
const scalarField& T,
const label patchi
) const;
//- Heat capacity at constant pressure/volume [J/kg/K]
virtual tmp<volScalarField> Cpv() const;
//- Heat capacity at constant pressure/volume for patch [J/kg/K]
virtual tmp<scalarField> Cpv
(
const scalarField& T,
const label patchi
) const;
// Fields derived from transport state variables
//- Thermal diffusivity for energy of mixture [kg/m/s]
virtual tmp<volScalarField> alphahe() const;
//- Thermal diffusivity for energy of mixture for patch [kg/m/s]
virtual tmp<scalarField> alphahe(const label patchi) const;
//- Effective thermal turbulent conductivity
// of mixture [W/m/K]
virtual tmp<volScalarField> kappaEff
(
const volScalarField&
) const;
//- Effective thermal turbulent diffusivity of mixture [kg/m/s]
virtual tmp<volScalarField> alphaEff
(
const volScalarField& alphat
) const;
//- Effective thermal turbulent conductivity
// of mixture for patch [W/m/K]
virtual tmp<scalarField> kappaEff
(
const scalarField& alphat,
const label patchi
) const;
//- Effective thermal turbulent diffusivity of mixture
// for patch [kg/m/s]
virtual tmp<scalarField> alphaEff
(
const scalarField& alphat,
const label patchi
) const;
//- Return true if thermal conductivity is isotropic
virtual bool isotropic() const
{
return true;
}
//- Anisotropic thermal conductivity [W/m/K]
virtual tmp<volVectorField> Kappa() const;
//- Anisotropic thermal conductivity [W/m/K]
virtual tmp<vectorField> Kappa(const label patchi) const;
//- Anisotropic thermal conductivity for patch
// in the local coordinate system [W/m/K]
virtual tmp<symmTensorField> KappaLocal(const label patchi) const;
//- Return the heat flux
virtual tmp<surfaceScalarField> q() const;
//- Return the source term for the energy equation
virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
//- Update properties
virtual void correct();
// I-O
//- Read thermophysicalProperties dictionary
virtual bool read();
}; };

View File

@ -37,7 +37,8 @@ boundaryField
} }
hole hole
{ {
type zeroGradient; type fixedValue;
value uniform 400;
} }
frontAndBack frontAndBack
{ {

View File

@ -32,7 +32,7 @@ E
value 2e+11; value 2e+11;
} }
Cp Cv
{ {
type uniform; type uniform;
value 434; value 434;

View File

@ -27,8 +27,6 @@ ddtSchemes
gradSchemes gradSchemes
{ {
default leastSquares; default leastSquares;
grad(D) leastSquares;
grad(T) leastSquares;
} }
divSchemes divSchemes
@ -39,9 +37,7 @@ divSchemes
laplacianSchemes laplacianSchemes
{ {
default none; default Gauss linear corrected;
laplacian(DD,D) Gauss linear corrected;
laplacian(kappa,T) Gauss linear corrected;
} }
interpolationSchemes interpolationSchemes

View File

@ -16,7 +16,7 @@ FoamFile
solvers solvers
{ {
"(D|T)" "(D|e)"
{ {
solver GAMG; solver GAMG;
tolerance 1e-06; tolerance 1e-06;

View File

@ -32,7 +32,7 @@ E
value 2e+11; value 2e+11;
} }
Cp Cv
{ {
type uniform; type uniform;
value 434; value 434;